CASM  1.1.0
A Clusters Approach to Statistical Mechanics
SupercellSymInfo.hh
Go to the documentation of this file.
1 #ifndef CASM_SupercellSymInfo
2 #define CASM_SupercellSymInfo
3 
4 #include <vector>
5 
10 #include "casm/global/eigen.hh"
14 
15 namespace CASM {
16 namespace SymRepTools {
17 struct IrrepInfo;
18 }
19 
20 namespace xtal {
21 class UnitCell;
22 }
23 using xtal::UnitCell;
24 
25 class PermuteIterator;
26 
32  public:
34  using SublatSymReps = std::vector<SymGroupRep::RemoteHandle>;
35 
38  SupercellSymInfo(Lattice const &_prim_lat, Lattice const &_super_lat,
39  Index number_of_sublats, SymGroup const &_prim_factor_group,
40  SymGroupRepID basis_permutation_symrep_ID,
41  std::map<DoFKey, SymGroupRepID> const &global_dof_symrep_IDs,
42  std::vector<SymGroupRepID> const &occ_symrep_IDs,
43  std::map<DoFKey, std::vector<SymGroupRepID> > const
44  &local_dof_symrep_IDs);
45 
49 
53 
56  std::map<DoFKey, SymGroupRep::RemoteHandle> const &global_dof_symreps()
57  const {
58  return m_global_dof_symreps;
59  }
60 
64  return m_global_dof_symreps.at(_key);
65  }
66 
76  SublatSymReps const &occ_symreps() const { return m_occ_symreps; }
77 
80  std::map<DoFKey, SublatSymReps> const &local_dof_symreps() const {
81  return m_local_dof_symreps;
82  }
83 
90  SublatSymReps const &local_dof_symreps(DoFKey const &_key) const {
91  return m_local_dof_symreps.at(_key);
92  }
93 
98  }
99 
103  const {
105  }
106 
111  const std::vector<Permutation> &translation_permutations() const {
113  }
114 
117  SymGroup const &factor_group() const { return m_factor_group; }
118 
121  bool has_aniso_occs() const { return m_has_aniso_occs; }
122 
125 
129  }
130 
132  const xtal::Lattice &prim_lattice() const {
134  }
135 
139  }
140 
147  }
148 
152 
155 
158  Index supercell_factor_group_index) const;
159 
160  // begin and end iterators for iterating over translation and factor group
161  // permutations
164  permute_const_iterator permute_it(Index supercell_factor_group_index,
165  Index translation_index) const;
166  permute_const_iterator permute_it(Index supercell_factor_group_index,
167  UnitCell translation) const;
168 
169  private:
173 
174  // TODO: I don't think this belongs in SupercellSymInfo, but neither did
175  // PrimGrid. I'm keeping the functionality where I found it for now, but we
176  // should consider moving it elsewhere
180 
181  // TODO: See TODO comment for m_unitcell_to_index_converter
185 
188  std::vector<Permutation> m_translation_permutations;
189 
190  // m_factor_group is factor group of the super cell, found by identifying the
191  // subgroup of
192  // (*this).prim().factor_group() that leaves the supercell lattice vectors
193  // unchanged if (*this).prim() is actually primitive, then
194  // m_factor_group.size() <= 48 NOTE: This is different from the SymGroup found
195  // by doing (*this).superstruc().factor_group()
196  // if Tprim is the translation group formed by the primitive cell
197  // lattice vectors, then m_factor_group is the group formed by the
198  // cosets of Tprim in the supercell space group if Tsuper is the
199  // translation group formed by the supercell lattice vectors, then,
200  // m_occupation(init_config.occupation()),
201  // m_displacement(init_config.displacement()),
202  // m_strain(init_config.supercell().strain
203  // (*this).superstruc().factor_group() is the group formed by the cosets
204  // of Tsuper in the supercell space group
206 
208 
210 
211  std::map<DoFKey, SublatSymReps> m_local_dof_symreps;
212 
213  std::map<DoFKey, SymGroupRep::RemoteHandle> m_global_dof_symreps;
214 
215  // true if there species with non-identity symreps
217 
218  // true if any site has occupation DoFs
220 
221  // m_perm_symrep_ID is the ID of the SymGroupRep of prim().factor_group() that
222  // describes how operations of m_factor_group permute sites of the Supercell.
223  // NOTE: The permutation representation is for (*this).prim().factor_group(),
224  // which may contain
225  // more operations than m_factor_group, so the Permutation SymGroupRep
226  // may have 'gaps' at the operations that aren't in m_factor_group. You
227  // should access elements of the SymGroupRep using the the
228  // Supercel::factor_group_permute(int) method, so that you don't
229  // encounter the gaps OR, see note for Supercell::permutation_symrep()
230  // below.
232 };
233 
234 std::string hermite_normal_form_name(const Eigen::Matrix3l &matrix);
235 
237 
239 std::string make_supercell_name(SymGroup const &point_group,
240  Lattice const &prim_lattice,
241  Lattice const &supercell_lattice);
242 
244 std::string make_canonical_supercell_name(SymGroup const &point_group,
245  Lattice const &prim_lattice,
246  Lattice const &supercell_lattice);
247 
250  SymGroup const &factor_group, Lattice const &prim_lattice,
251  std::string supercell_name);
252 
255 std::pair<MasterSymGroup, SymGroupRepID> make_collective_dof_symrep(
256  std::set<Index> const &site_indices, SupercellSymInfo const &_syminfo,
257  DoFKey const &_key, std::vector<PermuteIterator> const &_group);
258 
259 } // namespace CASM
260 
261 #endif
A class that collects all symmetry information for for performing symmetry transformations on the sit...
const xtal::Superlattice & superlattice() const
const reference to superlattice
SymGroupRep::RemoteHandle m_basis_perm_symrep
const Permutation & factor_group_permute(Index supercell_factor_group_index) const
Site permutation corresponding to supercell factor group operation.
const xtal::UnitCellCoordIndexConverter & unitcellcoord_index_converter() const
UnitCellCoordIndexConverter for this superstructure/primstructure pair Used to convert from lattice t...
SupercellSymInfo(Lattice const &_prim_lat, Lattice const &_super_lat, Index number_of_sublats, SymGroup const &_prim_factor_group, SymGroupRepID basis_permutation_symrep_ID, std::map< DoFKey, SymGroupRepID > const &global_dof_symrep_IDs, std::vector< SymGroupRepID > const &occ_symrep_IDs, std::map< DoFKey, std::vector< SymGroupRepID > > const &local_dof_symrep_IDs)
Construct with primitive and super lattice, number of sublattice and all relevant representation IDs.
xtal::UnitCellCoordIndexConverter m_unitcellcoord_to_index_converter
permute_const_iterator translate_end() const
End PermuteIterator over pure translational permutations.
const xtal::Lattice & prim_lattice() const
const reference to primitive lattice
SymGroup const & factor_group() const
Subgroup of primitive-cell factor group operations that leave supercell lattice invariant.
permute_const_iterator permute_end() const
const xtal::UnitCellIndexConverter & unitcell_index_converter() const
UnitCellIndexConverter for this superlattice/primlattice pair Used to convert from lattice translatio...
SymGroupRep::RemoteHandle const & site_permutation_symrep() const
const xtal::Lattice & supercell_lattice() const
const reference to supercell lattice
const std::vector< Permutation > & translation_permutations() const
Permutations describing reordering of sites of supercell due to a lattice translation of the primitiv...
SublatSymReps const & local_dof_symreps(DoFKey const &_key) const
SymGroupRep handle for site DoF matrix representation of supercell's factor group If sites are also p...
xtal::UnitCellIndexConverter m_unitcell_to_index_converter
PermuteIterator permute_const_iterator
bool has_occupation_dofs() const
true if any sublattice has more than one allowed occupant
std::vector< SymGroupRep::RemoteHandle > SublatSymReps
permute_const_iterator permute_begin() const
permute_const_iterator permute_it(Index supercell_factor_group_index, Index translation_index) const
std::map< DoFKey, SymGroupRep::RemoteHandle > m_global_dof_symreps
SymGroupRep::RemoteHandle const & basis_permutation_symrep() const
std::map< DoFKey, SymGroupRep::RemoteHandle > const & global_dof_symreps() const
Const reference global DoF matrix representations of the supercell's factor group.
permute_const_iterator translate_begin() const
Begin PermuteIterator over pure translational permutations / Equivalent to permute_begin()
Eigen::Matrix3l transformation_matrix_to_super() const
long-int transformation from primitive lattice vectors to supercell lattice vectors supercell_lattice...
std::map< DoFKey, SublatSymReps > m_local_dof_symreps
SymGroupRep::RemoteHandle m_site_perm_symrep
std::map< DoFKey, SublatSymReps > const & local_dof_symreps() const
Const reference local DoF matrix representations of the supercell's factor group.
bool has_aniso_occs() const
true if any species are anisotropic, such that the occ_symreps are non-trivial
xtal::Superlattice m_supercell_superlattice
std::vector< Permutation > m_translation_permutations
SublatSymReps const & occ_symreps() const
SymGroupRep handle for occupant permutation representation of supercell's factor group An occupant pe...
SymGroupRep::RemoteHandle const & global_dof_symrep(DoFKey const &_key) const
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:42
Type-safe ID object for communicating and accessing Symmetry representation info.
const Lattice & prim_lattice() const
Definition: Superlattice.hh:20
const Eigen::Matrix3l & transformation_matrix_to_super() const
Definition: Superlattice.hh:24
const Lattice & superlattice() const
Definition: Superlattice.hh:18
Unit Cell Indices.
xtal::Superlattice make_superlattice_from_supercell_name(Structure const &prim, std::string supercell_name)
Construct a Superlattice from the supercell name.
Definition: Supercell.cc:350
std::string make_canonical_supercell_name(Structure const &prim, xtal::Superlattice const &superlattice)
Make the canonical supercell name from a Superlattice.
Definition: Supercell.cc:343
std::string make_supercell_name(Structure const &prim, xtal::Superlattice const &superlattice)
Make the supercell name from a Superlattice.
Definition: Supercell.cc:336
SharedPrimFormatter< jsonParser > factor_group()
VectorXdSupercellSymInfoFormatter supercell_lattice()
Main CASM namespace.
Definition: APICommand.hh:8
std::pair< MasterSymGroup, SymGroupRepID > make_collective_dof_symrep(std::set< Index > const &site_indices, SupercellSymInfo const &_syminfo, DoFKey const &_key, std::vector< PermuteIterator > const &_group)
Make the matrix representation for group '_group' describing the transformation of DoF '_key' among a...
std::string DoFKey
Definition: DoFDecl.hh:7
Eigen::Matrix3l make_hermite_normal_form(std::string hermite_normal_form_name)
std::string hermite_normal_form_name(const Eigen::Matrix3l &matrix)
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Matrix< long int, 3, 3 > Matrix3l
Definition: eigen.hh:12