Hold information about a crystal structure.
More...
#include <structures.hpp>
|
| Crystal_structure (Eigen::Matrix3d _lattvec, Eigen::Matrix< double, 3, Eigen::Dynamic > _positions, std::vector< std::string > _elements, std::vector< int > _numbers) |
| Basic constructor.
|
|
| Crystal_structure (Eigen::Matrix3d _lattvec, Eigen::Matrix< double, 3, Eigen::Dynamic > _positions, std::vector< std::string > _elements, std::vector< int > _numbers, std::vector< double > _masses) |
| Specialized constructor allowing for custom masses.
|
|
std::vector< std::string >::size_type | get_nelements () const |
| Return the number of elements in the structure.
|
|
int | get_natoms () const |
| Return the number of atoms in the motif.
|
|
std::string | get_element (int i) const |
| Return the chemical symbol for an atom. More...
|
|
double | get_mass (int i) const |
| Return the mass of the i-th atom. More...
|
|
double | get_gfactor (int i) const |
| Return the g factor of the i-th atom. More...
|
|
bool | is_alloy () const |
| Check if any atom belongs to a virtual element. More...
|
|
Eigen::MatrixXd | map_to_firstbz (const Eigen::Ref< const Eigen::Vector3d > &q) const |
| Find all images of a q point in the first Brillouin zone. More...
|
|
|
const Eigen::Matrix3d | lattvec |
| Lattice vectors in nm.
|
|
const Eigen::Matrix< double, 3, Eigen::Dynamic > | positions |
| Atomic positions in lattice coordinates.
|
|
const std::vector< std::string > | elements |
| Elements present in the structure.
|
|
const double | V |
| Volume of the unit cell.
|
|
const Eigen::MatrixXd | rlattvec |
| Lattice vectors of the reciprocal lattice.
|
|
const std::vector< int > | numbers |
| Number of atoms of each element. More...
|
|
Hold information about a crystal structure.
◆ get_element()
std::string alma::Crystal_structure::get_element |
( |
int |
i | ) |
const |
|
inline |
Return the chemical symbol for an atom.
- Parameters
-
- Returns
- the chemical symbol
◆ get_gfactor()
double alma::Crystal_structure::get_gfactor |
( |
int |
i | ) |
const |
|
inline |
Return the g factor of the i-th atom.
- Parameters
-
- Returns
- the Pearson deviation coefficient of the i-th atom's mass
◆ get_mass()
double alma::Crystal_structure::get_mass |
( |
int |
i | ) |
const |
|
inline |
Return the mass of the i-th atom.
- Parameters
-
- Returns
- the mass in a.m.u
◆ is_alloy()
bool alma::Crystal_structure::is_alloy |
( |
| ) |
const |
|
inline |
Check if any atom belongs to a virtual element.
- Returns
- true is the structure describes an alloy, or false otherwise.
◆ map_to_firstbz()
Eigen::MatrixXd alma::Crystal_structure::map_to_firstbz |
( |
const Eigen::Ref< const Eigen::Vector3d > & |
q | ) |
const |
|
inline |
Find all images of a q point in the first Brillouin zone.
A point on the surface of the first BZ will have more than one possible image.
- Parameters
-
[in] | - | original q point in Cartesian coordinates |
- Returns
- the Cartesian coordinates of all images of the point in the first Brillouin zone, as columns of a matrix.
◆ numbers
const std::vector<int> alma::Crystal_structure::numbers |
Number of atoms of each element.
The first numbers[0] atoms belong to elements[0], the next numbers[1] atoms belong to elements[1], etc.
The documentation for this class was generated from the following file: