41 #include <boost/mpi.hpp> 42 #include <boost/math/distributions/normal.hpp> 43 #include <boost/serialization/array.hpp> 46 #include <boost/math/distributions/normal.hpp> 59 namespace serialization {
60 template <
class Archive>
63 const unsigned int file_version);
71 friend class boost::serialization::access;
72 template <
class Archive>
73 friend void boost::serialization::save_construct_data(
76 const unsigned int file_version);
79 const std::string& description,
83 const std::vector<Threeph_process>& processes,
84 const boost::mpi::communicator& comm);
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);
99 const double gaussian;
105 const std::array<std::size_t, 2>
q;
107 const std::array<std::size_t, 2>
alpha;
110 const std::array<std::size_t, 2>& _q,
111 const std::array<std::size_t, 2>& _alpha,
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)) {
122 : domega(original.domega), sigma(original.sigma),
123 gaussian(original.gaussian), c(original.c), q(original.q),
124 alpha(original.alpha) {
132 const std::array<std::size_t, 2>& _q,
133 const std::array<std::size_t, 2>& _alpha,
137 : domega(_domega), sigma(_sigma), gaussian(_gaussian), c(_c),
138 q(
std::move(_q)), alpha(
std::move(_alpha)) {
159 const Eigen::Ref<const Eigen::VectorXd>& gfactors,
174 const boost::mpi::communicator& communicator,
175 double scalebroad = 1.0);
187 const std::vector<alma::Twoph_process>& processes,
188 const boost::mpi::communicator& comm);
202 const Eigen::Ref<const Eigen::VectorXd>& gfactors,
204 const std::vector<alma::Twoph_process>& processes,
205 const boost::mpi::communicator& comm);
209 namespace serialization {
214 template <
class Archive>
217 const unsigned int file_version) {
219 ar << make_array(t->
q.data(), t->
q.size());
220 ar << make_array(t->
alpha.data(), t->
alpha.size());
232 template <
class Archive>
235 const unsigned int file_version) {
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());
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
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