CASM  1.1.0
A Clusters Approach to Statistical Mechanics
SimpleStructure.hh
Go to the documentation of this file.
1 #ifndef SIMPLESTRUCTURE_HH
2 #define SIMPLESTRUCTURE_HH
3 
4 #include <map>
5 #include <set>
6 #include <string>
7 #include <vector>
8 
9 #include "casm/external/Eigen/Dense"
11 
12 namespace CASM {
13 namespace xtal {
14 
19 namespace SimpleStructureTools {
20 // TODO: Capitalize so it looks the same as COORD_MODE? Put it in some
21 // crystallography declarations file?
24 enum class SpeciesMode { ATOM, MOL };
25 } // namespace SimpleStructureTools
26 
33  public:
35 
44  struct Info {
47  using Coord = Eigen::MatrixXd::ColXpr;
48  using ConstCoord = Eigen::MatrixXd::ConstColXpr;
49 
51  std::vector<std::string> names;
52 
56 
60  std::map<std::string, Eigen::MatrixXd> properties;
61 
65  std::vector<Index> sort_by_name();
66 
68  Coord cart_coord(Index i) { return coords.col(i); }
69 
71  ConstCoord cart_coord(Index i) const { return coords.col(i); }
72 
75  void resize(Index N) {
76  names.resize(N, "Va");
77  coords.setZero(3, N);
78  for (auto &el : properties) el.second.setZero(Eigen::NoChange, N);
79  }
80 
82  Index size() const { return names.size(); }
83  };
84 
85  Index n_mol() const { return mol_info.size(); }
86 
87  Index n_atom() const { return atom_info.size(); }
88 
90  Info &info(SpeciesMode _mode) {
91  return _mode == SpeciesMode::MOL ? mol_info : atom_info;
92  }
93 
95  Info const &info(SpeciesMode _mode) const {
96  return _mode == SpeciesMode::MOL ? mol_info : atom_info;
97  }
98 
100  void within(double tol = TOL);
101 
104  void deform_coords(Eigen::Ref<const Eigen::Matrix3d> const &_F);
105 
108  void rotate_coords(Eigen::Ref<const Eigen::Matrix3d> const &_R);
109 
111 
112  // Use occupation vector in order to avoid messy molecule-name aliasing issues
114 
116 
117  std::map<std::string, Eigen::MatrixXd> properties;
118 };
119 
121 } // namespace xtal
122 } // namespace CASM
123 
124 #endif
Representation of a crystal of molecular and/or atomic occupants, and any additional properties....
void deform_coords(Eigen::Ref< const Eigen::Matrix3d > const &_F)
Apply homogeneous deformation gradient tensor _F to lat_column_mat, mol_info, and atom_info.
Info const & info(SpeciesMode _mode) const
Const info of specified context (atomic/molecular)
std::map< std::string, Eigen::MatrixXd > properties
void within(double tol=TOL)
Map all coordinates within the unit cell.
void rotate_coords(Eigen::Ref< const Eigen::Matrix3d > const &_R)
Apply homogeneous rotation matrix _R to lat_column_mat, mol_info, and atom_info.
Info & info(SpeciesMode _mode)
Info of specified context (atomic/molecular)
SpeciesMode
enum to refer to a particular representation of the occupants (atomic or molecular)
Main CASM namespace.
Definition: APICommand.hh:8
Eigen::MatrixXd MatrixXd
Eigen::Matrix3d Matrix3d
const double TOL
Definition: definitions.hh:30
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Struct to encode all information about the crystal basis Info may describe the basis in a atomic cont...
Eigen::MatrixXd coords
(3 x names.size()) matrix of coordinates. coords.col(i) is Cartesian coordinate of site 'i'
Eigen::MatrixXd::ColXpr Coord
ConstCoord cart_coord(Index i) const
const access for coordinate of site 'i'
std::vector< std::string > names
names[i] is name of species that occupies sites 'i'
Coord cart_coord(Index i)
Access coordinate of site 'i'.
Index size() const
Number of sites is defined as names.size()
std::map< std::string, Eigen::MatrixXd > properties
map of [property name, (m x names.size()) matrix] for all numerical site properties properties are as...
std::vector< Index > sort_by_name()
permutation that results in sites sorted alphabetically by species Guaranteed stable: will not change...
void resize(Index N)
Resize to hold N sites. All coordinates set to zero, all occupants set to "Va".
Eigen::MatrixXd::ConstColXpr ConstCoord