CASM  1.1.0
A Clusters Approach to Statistical Mechanics
Supercell.hh
Go to the documentation of this file.
1 #ifndef CASM_Supercell
2 #define CASM_Supercell
3 
8 #include "casm/database/Named.hh"
14 
15 namespace CASM {
16 
17 namespace xtal {
18 class Site;
19 class Coordinate;
20 class UnitCellCoord;
21 class IntegralCoordinateWithin_f;
22 } // namespace xtal
23 using xtal::Coordinate;
24 using xtal::Site;
25 using xtal::UnitCellCoord;
26 
27 template <typename T, typename U>
29 class PermuteIterator;
30 class PrimClex;
31 class Clexulator;
32 class PrimNeighborList;
33 class SuperNeighborList;
34 class Structure;
35 
36 namespace DB {
37 template <typename T>
38 class DatabaseIterator;
39 }
40 
49 class Supercell
50  : public DB::Named<
51  Comparisons<SupercellCanonicalForm<CRTPBase<Supercell>>>> {
52  public:
54 
55  // **** Constructors ****
56 
57  Supercell(const Supercell &RHS);
58 
59  Supercell(std::shared_ptr<Structure const> const &_shared_prim,
60  const Lattice &superlattice);
61  Supercell(std::shared_ptr<Structure const> const &_shared_prim,
62  Eigen::Matrix3l const &superlattice_matrix);
63 
64  Supercell(const PrimClex *_prim, const Lattice &superlattice);
65  Supercell(const PrimClex *_prim,
66  const Eigen::Ref<const Eigen::Matrix3l> &superlattice_matrix);
67 
68  ~Supercell();
69 
70  // **** Coordinates ****
71 
74 
77  Index linear_index(const Coordinate &coord, double tol = TOL) const;
78 
80  Index linear_index(const UnitCellCoord &bijk) const;
81 
84 
87 
90  Eigen::VectorXi max_allowed_occupation() const;
91 
92  // **** Accessors ****
93 
94  const Structure &prim() const;
95 
96  std::shared_ptr<Structure const> const &shared_prim() const;
97 
98  double crystallography_tol() const;
99 
101  bool has_primclex() const;
102 
104  void set_primclex(PrimClex const *primclex_ptr) const;
105 
107  const PrimClex &primclex() const;
108 
110  Index volume() const;
111 
112  Index basis_size() const;
113 
114  Index num_sites() const;
115 
117 
118  Eigen::Matrix3l transf_mat() const;
119 
121  const Lattice &lattice() const;
122 
134  const SuperNeighborList &nlist() const;
135 
136  // Factor group of this supercell
137  const SymGroup &factor_group() const;
138 
139  // SymInfo object of this supercell
140  const SupercellSymInfo &sym_info() const;
141 
142  bool operator<(const Supercell &B) const;
143 
145  std::pair<DB::DatabaseIterator<Supercell>, bool> insert() const;
146 
147  private:
150 
151  bool eq_impl(const Supercell &B) const;
152 
153  // **** Generating functions ****
154 
155  std::string generate_name_impl() const;
156 
157  // Note:
158  // - Prefer not to access PrimClex via Supercell. In future, PrimClex access
159  // via Supercell will
160  // be removed completely.
161  // - Until this is removed, it may be nullptr, in which case, some features
162  // will throw. Only
163  // access via this->primclex() so that an error will be thrown if m_primclex
164  // is nullptr.
165  // - Mutable as a temporary workaround
166  mutable PrimClex const *m_primclex;
167 
168  std::shared_ptr<Structure const> m_shared_prim;
169 
171 
174 
179 };
180 
182 std::string make_supercell_name(Structure const &prim,
183  xtal::Superlattice const &superlattice);
184 
187  Structure const &prim, xtal::Superlattice const &superlattice);
188 
191  Structure const &prim, std::string supercell_name);
192 
194 Supercell &apply(const SymOp &op, Supercell &scel);
195 
197 Supercell copy_apply(const SymOp &op, const Supercell &scel);
198 
199 // --- The following are deprecated ----
200 
201 const Supercell &make_supercell(const PrimClex &primclex, std::string name);
202 
203 std::shared_ptr<Supercell> make_shared_supercell(const PrimClex &primclex,
204  std::string name);
205 
206 Eigen::Matrix3l transf_mat(const Lattice &prim_lat, const Lattice &super_lat);
207 
208 std::string generate_name(const Eigen::Matrix3l &transf_mat);
209 
210 std::string scelname(const Structure &prim, const Lattice &superlat);
211 
212 std::string canonical_scelname(const Structure &prim, const Lattice &superlat);
213 
214 namespace xtal {
215 IntegralCoordinateWithin_f make_bring_within_f(const Supercell &scel);
216 }
217 
219 } // namespace CASM
220 #endif
Evaluates correlations.
Definition: Clexulator.hh:440
PrimClex is the top-level data structure for a CASM project.
Definition: PrimClex.hh:55
The PrimNeighborList gives the coordinates of UnitCell that are neighbors of the origin UnitCell.
Definition: NeighborList.hh:39
Structure specifies the lattice and atomic basis of a crystal.
Definition: Structure.hh:30
Represents a supercell of the primitive parent crystal structure.
Definition: Supercell.hh:51
Index sublat(Index linear_index) const
Return the sublattice index for a linear index.
Definition: Supercell.cc:165
Index linear_index(const Coordinate &coord, double tol=TOL) const
Given a Coordinate and tolerance, return linear index into Configuration.
Definition: Supercell.cc:183
SupercellSymInfo m_sym_info
Definition: Supercell.hh:170
UnitCellCoord uccoord(Index linear_index) const
Return the integral coordinates corresponding to a linear index.
Definition: Supercell.cc:209
Eigen::VectorXi max_allowed_occupation() const
returns maximum allowed occupation bitstring – used for initializing enumeration counters
Definition: Supercell.cc:213
SymGroupRep const & permutation_symrep() const
const PrimClex & primclex() const
Use while transitioning Supercell to no longer need a PrimClex const *
Definition: Supercell.cc:150
double crystallography_tol() const
Definition: Supercell.cc:119
Coordinate coord(Index linear_index) const
Return the linear index corresponding to integral coordinates.
Definition: Supercell.cc:199
bool has_primclex() const
Use while transitioning Supercell to no longer need a PrimClex const *
Definition: Supercell.cc:127
bool operator<(const Supercell &B) const
Definition: Supercell.cc:267
Supercell(const Supercell &RHS)
Definition: Supercell.cc:36
const Lattice & lattice() const
The super lattice.
Definition: Supercell.cc:239
const SuperNeighborList & nlist() const
Returns the SuperNeighborList.
Definition: Supercell.cc:244
void set_primclex(PrimClex const *primclex_ptr) const
Use while transitioning Supercell to no longer need a PrimClex const *
Definition: Supercell.cc:136
Eigen::Matrix3l transf_mat() const
Definition: Supercell.cc:235
Index basis_size() const
Definition: Supercell.cc:231
std::pair< DB::DatabaseIterator< Supercell >, bool > insert() const
Insert the canonical form of this into the database.
Definition: Supercell.cc:283
PrimClex const * m_primclex
Definition: Supercell.hh:166
Index m_nlist_size_at_construction
Definition: Supercell.hh:178
std::shared_ptr< Structure const > m_shared_prim
Definition: Supercell.hh:168
notstd::cloneable_ptr< SuperNeighborList > m_nlist
SuperNeighborList, mutable for lazy construction.
Definition: Supercell.hh:173
Index num_sites() const
Definition: Supercell.cc:233
Index volume() const
Return number of primitive cells that fit inside of *this.
Definition: Supercell.cc:227
bool eq_impl(const Supercell &B) const
Definition: Supercell.cc:289
const SymGroup & factor_group() const
Definition: Supercell.cc:260
const Structure & prim() const
Definition: Supercell.cc:113
std::shared_ptr< Structure const > const & shared_prim() const
Definition: Supercell.cc:115
const SupercellSymInfo & sym_info() const
Definition: Supercell.cc:265
std::string generate_name_impl() const
Return supercell name.
Definition: Supercell.cc:317
A class that collects all symmetry information for for performing symmetry transformations on the sit...
PermuteIterator permute_const_iterator
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:42
SymGroupRep is an alternative representation of a SymGroup for something other than real space....
Definition: SymGroupRep.hh:31
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
Definition: SymOp.hh:28
Represents cartesian and fractional coordinates.
Definition: Coordinate.hh:34
Unit Cell Coordinates.
A 'cloneable_ptr' can be used in place of 'unique_ptr'.
std::shared_ptr< Supercell > make_shared_supercell(const PrimClex &primclex, std::string name)
Definition: Supercell.cc:433
std::string generate_name(const Eigen::Matrix3l &transf_mat)
Make hermite normal form name [deprecated].
Definition: Supercell.cc:482
const Supercell & make_supercell(const PrimClex &primclex, std::string name)
Definition: Supercell.cc:378
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 canonical_scelname(const Structure &prim, const Lattice &superlat)
Make canonical supercell name name [deprecated].
Definition: Supercell.cc:510
std::string scelname(const Structure &prim, const Lattice &superlat)
Make supercell name name [deprecated].
Definition: Supercell.cc:497
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
Supercell & apply(const SymOp &op, Supercell &scel)
Apply symmetry operation to Supercell.
Definition: Supercell.cc:357
Eigen::Matrix3l transf_mat(const Lattice &prim_lat, const Lattice &super_lat)
IntegralCoordinateWithin_f make_bring_within_f(const Supercell &scel)
Make IntegralCoordinateWithin_f for Supercell [deprecated].
Definition: Supercell.cc:520
Main CASM namespace.
Definition: APICommand.hh:8
MappingNode copy_apply(PermuteIterator const &_it, MappingNode const &_node, bool transform_cost_mat=true)
Reorders the permutation and compounds the spatial isometry (rotation.
const double TOL
Definition: definitions.hh:30
GenericDatumFormatter< std::string, DataObject > name()
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
Matrix< long int, 3, 3 > Matrix3l
Definition: eigen.hh:12
PrimClex * primclex
Definition: settings.cc:135
Implements other comparisons in terms of '<'.
Definition: Comparisons.hh:25