AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
dynamical_matrix.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 
19 #include <boost/serialization/complex.hpp>
20 #include <boost/serialization/array.hpp>
21 #include <structures.hpp>
22 #include <symmetry.hpp>
23 
24 namespace alma {
31 private:
32  friend class boost::serialization::access;
39  template <class Archive>
40  void save(Archive& ar, const unsigned int version) const {
41  auto rows = this->omega.size();
42 
43  ar << rows;
44  ar << boost::serialization::make_array(this->omega.data(), rows);
45  rows = this->wfs.rows();
46  auto cols = this->wfs.cols();
47  ar << rows;
48  ar << cols;
49  ar << boost::serialization::make_array(this->wfs.data(), rows * cols);
50  rows = this->vg.rows();
51  cols = this->vg.cols();
52  ar << rows;
53  ar << cols;
54  ar << boost::serialization::make_array(this->vg.data(), rows * cols);
55  }
56 
57 
64  template <class Archive>
65  void load(Archive& ar, const unsigned int version) {
66  decltype(this->omega.size()) rows;
67  ar >> rows;
68  this->omega.resize(rows);
69  ar >> boost::serialization::make_array(this->omega.data(), rows);
70  decltype(rows) cols;
71  ar >> rows;
72  ar >> cols;
73  this->wfs.resize(rows, cols);
74  ar >> boost::serialization::make_array(this->wfs.data(), rows * cols);
75  ar >> rows;
76  ar >> cols;
77  this->vg.resize(rows, cols);
78  ar >> boost::serialization::make_array(this->vg.data(), rows * cols);
79  }
80 
81 
82  // Instruct boost::serialization to expect separate load()
83  // and save() methods.
84  BOOST_SERIALIZATION_SPLIT_MEMBER()
85 public:
88  Eigen::ArrayXd omega;
91  Eigen::MatrixXcd wfs;
93  Eigen::ArrayXXd vg;
95  Spectrum_at_point(const Eigen::Ref<const Eigen::ArrayXd>& _omega,
96  const Eigen::Ref<const Eigen::MatrixXcd>& _wfs,
97  const Eigen::Ref<const Eigen::ArrayXXd>& _vg)
98  : omega(_omega), wfs(_wfs), vg(_vg) {
99  }
100 
101 
104  }
105 };
106 
109 public:
112  const Symmetry_operations& syms,
113  const Harmonic_ifcs& fcs);
116  const Symmetry_operations& syms,
117  const Harmonic_ifcs& fcs,
118  const Dielectric_parameters& born);
125  std::array<Eigen::MatrixXcd, 4> build(
126  const Eigen::Ref<const Eigen::Vector3d>& q) const;
127 
134  std::unique_ptr<Spectrum_at_point> get_spectrum(
135  const Eigen::Ref<const Eigen::Vector3d>& q) const;
136 
142  std::size_t get_nblocks() const {
143  return this->blocks.size();
144  }
150  Eigen::MatrixXd get_block(std::size_t i) const {
151  if (i >= this->blocks.size())
152  throw value_error("invalid index");
153  return this->blocks[i];
154  }
155 
156 
162  Triple_int get_pos(std::size_t i) const {
163  if (i >= this->blocks.size())
164  throw value_error("invalid index");
165  return this->pos[i];
166  }
167 
168 
169 private:
173  std::vector<Eigen::MatrixXd> blocks;
176  std::vector<Eigen::MatrixXd> masks;
178  std::vector<Triple_int> pos;
181  const int na;
184  const int nb;
187  const int nc;
189  const double V;
191  const Crystal_structure structure;
193  const Eigen::MatrixXd rlattvec;
195  const Symmetry_operations symmetries;
197  Eigen::MatrixXd mpos;
199  Eigen::MatrixXd cpos;
203  Eigen::ArrayXXd massmatrix;
206  const bool nonanalytic;
209  const Dielectric_parameters born;
216  void copy_blocks(const Harmonic_ifcs& fcs);
217 
223  Eigen::ArrayXcd get_exponentials(
224  const Eigen::Ref<const Eigen::Vector3d>& q) const;
225 
233  std::array<Eigen::ArrayXXd, 4> build_nac(
234  const Eigen::Ref<const Eigen::Vector3d>& q) const;
235 };
236 } // namespace alma
std::size_t get_nblocks() const
Obtain the number of building blocks of unfolded, mass-reduced interatomic force constants used to bu...
Definition: dynamical_matrix.hpp:142
Definition: analytic1d.hpp:26
Spectrum_at_point(const Eigen::Ref< const Eigen::ArrayXd > &_omega, const Eigen::Ref< const Eigen::MatrixXcd > &_wfs, const Eigen::Ref< const Eigen::ArrayXXd > &_vg)
Basic constructor.
Definition: dynamical_matrix.hpp:95
Spectrum_at_point()
Empty default constructor.
Definition: dynamical_matrix.hpp:103
Exception related to the parameters passed to a function.
Definition: exceptions.hpp:33
Hold information about the harmonic interactions between atoms.
Definition: structures.hpp:233
std::array< int, 3 > Triple_int
Convenient shorthand for an array of three ints.
Definition: structures.hpp:43
Eigen::ArrayXXd vg
Cartesian components of the group velocities in km / s.
Definition: dynamical_matrix.hpp:93
Hold information about the polarization properties of the structure.
Definition: structures.hpp:272
Eigen::ArrayXd omega
Angular frequencies, in rad / ps.
Definition: dynamical_matrix.hpp:88
Factory of Dynamical_matrix objects.
Definition: dynamical_matrix.hpp:108
Eigen::MatrixXd get_block(std::size_t i) const
Obtain a copy of one of the building blocks of the dynamical matrix.
Definition: dynamical_matrix.hpp:150
Triple_int get_pos(std::size_t i) const
Obtain the direct coordinates of the lattice point associated to a block.
Definition: dynamical_matrix.hpp:162
Definitions of the basic data-handling classes in ALMA.
Objects of this class hold a subset of the information provided by spg_get_dataset().
Definition: symmetry.hpp:47
C++ interface to Atsushi Togo&#39;s spglib.
Hold information about a crystal structure.
Definition: structures.hpp:51
POD class that holds all the information about the harmonic properties of the system at a particular ...
Definition: dynamical_matrix.hpp:30
Eigen::MatrixXcd wfs
Eigenvectors, normalized on a single unit cell.
Definition: dynamical_matrix.hpp:91