AlmaBTE  1.3
A solver of the space- and time-dependent Boltzmann transport equation for phonons
isotopic_scattering.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 
36 
37 #include <cstddef>
38 #include <cmath>
39 #include <array>
40 #include <vector>
41 #include <boost/mpi.hpp>
42 #include <boost/math/distributions/normal.hpp>
43 #include <boost/serialization/array.hpp>
44 #include <structures.hpp>
45 #include <qpoint_grid.hpp>
46 #include <boost/math/distributions/normal.hpp>
47 
48 // Forward declarations of elements documented later on.
49 // Note that serialization of std::vectors of objects without
50 // a default constructor is broken in Boost version 1.58.0. See:
51 // http://stackoverflow.com/a/30437359/85371
52 // https://svn.boost.org/trac/boost/ticket/11342
53 // This will prevent ALMA from compiling. We check for that
54 // specific version in our cmake configuration.
55 namespace alma {
56 class Twoph_process;
57 }
58 namespace boost {
59 namespace serialization {
60 template <class Archive>
61 inline void save_construct_data(Archive& ar,
62  const alma::Twoph_process* t,
63  const unsigned int file_version);
64 }
65 } // namespace boost
66 
67 namespace alma {
70 private:
71  friend class boost::serialization::access;
72  template <class Archive>
73  friend void boost::serialization::save_construct_data(
74  Archive& ar,
75  const Twoph_process* t,
76  const unsigned int file_version);
77 
78  friend void save_bulk_hdf5(const char* filename,
79  const std::string& description,
80  const Crystal_structure& cell,
81  const Symmetry_operations& symmetries,
82  const Gamma_grid& grid,
83  const std::vector<Threeph_process>& processes,
84  const boost::mpi::communicator& comm);
85 
86  friend std::tuple<std::string,
87  std::unique_ptr<Crystal_structure>,
88  std::unique_ptr<Symmetry_operations>,
89  std::unique_ptr<Gamma_grid>,
90  std::unique_ptr<std::vector<Threeph_process>>>
91  load_bulk_hdf5(const char* filename, const boost::mpi::communicator& comm);
92 
94  const double domega;
96  const double sigma;
99  const double gaussian;
100 
101 public:
103  const std::size_t c;
105  const std::array<std::size_t, 2> q;
107  const std::array<std::size_t, 2> alpha;
109  Twoph_process(std::size_t _c,
110  const std::array<std::size_t, 2>& _q,
111  const std::array<std::size_t, 2>& _alpha,
112  double _domega,
113  double _sigma)
114  : domega(_domega), sigma(_sigma),
115  gaussian(boost::math::pdf(boost::math::normal(0., sigma), domega)),
116  c(_c), q(std::move(_q)), alpha(std::move(_alpha)) {
117  }
118 
119 
121  Twoph_process(const Twoph_process& original)
122  : domega(original.domega), sigma(original.sigma),
123  gaussian(original.gaussian), c(original.c), q(original.q),
124  alpha(original.alpha) {
125  }
126 
127 
131  Twoph_process(std::size_t _c,
132  const std::array<std::size_t, 2>& _q,
133  const std::array<std::size_t, 2>& _alpha,
134  double _domega,
135  double _sigma,
136  double _gaussian)
137  : domega(_domega), sigma(_sigma), gaussian(_gaussian), c(_c),
138  q(std::move(_q)), alpha(std::move(_alpha)) {
139  }
140 
141 
147  double compute_gamma(const Crystal_structure& cell,
148  const Gamma_grid& grid) const;
149 
158  double compute_gamma(const Crystal_structure& cell,
159  const Eigen::Ref<const Eigen::VectorXd>& gfactors,
160  const Gamma_grid& grid) const;
161 };
172 std::vector<Twoph_process> find_allowed_twoph(
173  const Gamma_grid& grid,
174  const boost::mpi::communicator& communicator,
175  double scalebroad = 1.0);
176 
185 Eigen::ArrayXXd calc_w0_twoph(const Crystal_structure& cell,
186  const alma::Gamma_grid& grid,
187  const std::vector<alma::Twoph_process>& processes,
188  const boost::mpi::communicator& comm);
189 
201 Eigen::ArrayXXd calc_w0_twoph(const Crystal_structure& cell,
202  const Eigen::Ref<const Eigen::VectorXd>& gfactors,
203  const alma::Gamma_grid& grid,
204  const std::vector<alma::Twoph_process>& processes,
205  const boost::mpi::communicator& comm);
206 } // namespace alma
207 
208 namespace boost {
209 namespace serialization {
214 template <class Archive>
215 inline void save_construct_data(Archive& ar,
216  const alma::Twoph_process* t,
217  const unsigned int file_version) {
218  ar << t->c;
219  ar << make_array(t->q.data(), t->q.size());
220  ar << make_array(t->alpha.data(), t->alpha.size());
221  ar << t->domega;
222  ar << t->sigma;
223  ar << t->gaussian;
224 }
225 
226 
232 template <class Archive>
233 inline void load_construct_data(Archive& ar,
235  const unsigned int file_version) {
236  std::size_t c;
237 
238  ar >> c;
239  std::array<std::size_t, 2> q;
240  ar >> make_array(q.data(), q.size());
241  std::array<std::size_t, 2> alpha;
242  ar >> make_array(alpha.data(), alpha.size());
243  double domega;
244  ar >> domega;
245  double sigma;
246  ar >> sigma;
247  double gaussian;
248  ar >> gaussian;
249  ::new (t) alma::Twoph_process(c, q, alpha, domega, sigma, gaussian);
250 }
251 } // namespace serialization
252 } // namespace boost
Twoph_process(std::size_t _c, const std::array< std::size_t, 2 > &_q, const std::array< std::size_t, 2 > &_alpha, double _domega, double _sigma, double _gaussian)
Explicit constructor used when deserializing objects of this class.
Definition: isotopic_scattering.hpp:131
Definition: analytic1d.hpp:26
Definition: isotopic_scattering.hpp:58
void save_construct_data(Archive &ar, const alma::Twoph_process *t, const unsigned int file_version)
Overload required to serialize the const members of alma::Twoph_process.
Definition: isotopic_scattering.hpp:215
const std::array< std::size_t, 2 > q
q point indices of each of the two phonons involved.
Definition: isotopic_scattering.hpp:105
Representation of a elastic two-phonon process.
Definition: isotopic_scattering.hpp:69
void save_bulk_hdf5(const char *filename, const std::string &description, const Crystal_structure &cell, const Symmetry_operations &symmetries, const Gamma_grid &grid, const std::vector< Threeph_process > &processes, const boost::mpi::communicator &comm)
Save all relevant information about a bulk material to an HDF5 file.
Definition: bulk_hdf5.cpp:125
STL namespace.
const std::size_t c
Equivalence class of the first phonon.
Definition: isotopic_scattering.hpp:103
const std::array< std::size_t, 2 > alpha
Mode indices of the two phonons involved.
Definition: isotopic_scattering.hpp:107
Classes and functions used to manipulate grids in reciprocal space.
Twoph_process(const Twoph_process &original)
Copy constructor.
Definition: isotopic_scattering.hpp:121
std::vector< Twoph_process > find_allowed_twoph(const Gamma_grid &grid, const boost::mpi::communicator &communicator, double scalebroad=1.0)
Look for allowed two-phonon processes in a regular grid.
Definition: isotopic_scattering.cpp:53
void load_construct_data(Archive &ar, alma::Twoph_process *t, const unsigned int file_version)
Overload required to deserialize the const members of alma::Twoph_process by calling the non-default ...
Definition: isotopic_scattering.hpp:233
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
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_w0_twoph(const Crystal_structure &cell, const alma::Gamma_grid &grid, const std::vector< alma::Twoph_process > &processes, const boost::mpi::communicator &comm)
Compute and the two-phonon contribution to the RTA scattering rates for all vibrational modes on a gr...
Definition: isotopic_scattering.cpp:158
Twoph_process(std::size_t _c, const std::array< std::size_t, 2 > &_q, const std::array< std::size_t, 2 > &_alpha, double _domega, double _sigma)
Basic constructor.
Definition: isotopic_scattering.hpp:109
std::tuple< std::string, std::unique_ptr< Crystal_structure >, std::unique_ptr< Symmetry_operations >, std::unique_ptr< Gamma_grid >, std::unique_ptr< std::vector< Threeph_process > > > load_bulk_hdf5(const char *filename, const boost::mpi::communicator &comm)
Reconstruct all data structures from a file saved with save_bulk_hdf5().
Definition: bulk_hdf5.cpp:442