19 #include <boost/serialization/complex.hpp> 20 #include <boost/serialization/array.hpp> 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();
44 ar << boost::serialization::make_array(this->
omega.data(), rows);
45 rows = this->
wfs.rows();
46 auto cols = this->
wfs.cols();
49 ar << boost::serialization::make_array(this->
wfs.data(), rows * cols);
50 rows = this->
vg.rows();
51 cols = this->
vg.cols();
54 ar << boost::serialization::make_array(this->
vg.data(), rows * cols);
64 template <
class Archive>
65 void load(Archive& ar,
const unsigned int version) {
66 decltype(this->
omega.size()) rows;
68 this->
omega.resize(rows);
69 ar >> boost::serialization::make_array(this->
omega.data(), rows);
73 this->
wfs.resize(rows, cols);
74 ar >> boost::serialization::make_array(this->
wfs.data(), rows * cols);
77 this->
vg.resize(rows, cols);
78 ar >> boost::serialization::make_array(this->
vg.data(), rows * cols);
84 BOOST_SERIALIZATION_SPLIT_MEMBER()
96 const Eigen::Ref<const Eigen::MatrixXcd>& _wfs,
97 const Eigen::Ref<const Eigen::ArrayXXd>& _vg)
98 : omega(_omega), wfs(_wfs), vg(_vg) {
125 std::array<Eigen::MatrixXcd, 4> build(
126 const Eigen::Ref<const Eigen::Vector3d>& q)
const;
134 std::unique_ptr<Spectrum_at_point> get_spectrum(
135 const Eigen::Ref<const Eigen::Vector3d>& q)
const;
143 return this->blocks.size();
151 if (i >= this->blocks.size())
153 return this->blocks[i];
163 if (i >= this->blocks.size())
173 std::vector<Eigen::MatrixXd> blocks;
176 std::vector<Eigen::MatrixXd> masks;
178 std::vector<Triple_int> pos;
193 const Eigen::MatrixXd rlattvec;
197 Eigen::MatrixXd mpos;
199 Eigen::MatrixXd cpos;
203 Eigen::ArrayXXd massmatrix;
206 const bool nonanalytic;
223 Eigen::ArrayXcd get_exponentials(
224 const Eigen::Ref<const Eigen::Vector3d>& q)
const;
233 std::array<Eigen::ArrayXXd, 4> build_nac(
234 const Eigen::Ref<const Eigen::Vector3d>& q)
const;
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'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