CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
BasicStructure.hh
Go to the documentation of this file.
1 #ifndef BASICSTRUCTURE_HH
2 #define BASICSTRUCTURE_HH
3 
4 #include <iostream>
5 #include <cstdlib>
6 #include <cmath>
7 
9 
10 namespace CASM {
11  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
12  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
13 
14  class SymGroupRepID;
15  class Coordinate;
16  class UnitCellCoord;
17  class SiteCluster;
18  class MasterSymGroup;
19  template<typename ClustType> class GenericOrbitree;
20  typedef GenericOrbitree<SiteCluster> SiteOrbitree;
21 
28  template<typename CoordType>
30  class BasicStructure {
31  protected:
33 
35  bool SD_flag;
36 
37  public: // PUBLIC DATA MEMBERS -- (long-term, at least lattice should be made private and only updated via Structure::set_lattice)
39  std::string title;
40 
43 
44 
45  private: // PRIVATE METHODS
46 
47  void main_print(std::ostream &stream, COORD_TYPE mode, bool version5, int option) const;
48 
49  public: // PUBLIC METHODS
50 
51  // ****Constructors****
52  BasicStructure(const Lattice &init_lat) : m_lattice(init_lat) {};
53  BasicStructure() : m_lattice() {}; //added by Ivy (do we need/want this??)
54  BasicStructure(const fs::path &filepath);
55 
58  BasicStructure(const BasicStructure &RHS);
59 
60  virtual ~BasicStructure() {}
61 
62  // ****Inspectors/Accessors****
63 
66  template<typename CoordType2>
67  Index find(const CoordType2 &test_site, double tol = TOL) const;
68 
71  template<typename CoordType2>
72  Index find(const CoordType2 &test_site, const Coordinate &shift, double tol) const;
73 
74  const Lattice &lattice() const {
75  return m_lattice;
76  }
77 
80  template<typename CoordType2>
81  UnitCellCoord get_unit_cell_coord(const CoordType2 &test_site, double tol = TOL)const;
82 
83  // ****Mutators****
84 
85  // - Basic assignment/bookkeeping
86 
89  virtual BasicStructure &operator=(const BasicStructure &RHS);
90 
94  void copy_attributes_from(const BasicStructure &RHS);
95  //virtual void copy_attributes_from(const BasicStructure &RHS); <---- should be virtual, but will have to use some sort of visitor pattern to make work
96 
97  // update does reset() first, and then calls set_site_internals
98  void update();
99 
100  // clears site_internals and does within()
101  virtual void reset();
102 
104  void set_site_internals();
105 
107  void within();
108 
109  CoordType get_site(const UnitCellCoord &ucc) const;
110 
114  void set_lattice(const Lattice &lattice, COORD_TYPE mode);
115 
117  void set_basis(Array<CoordType> basis_in);
118 
119  // - Symmetry
120 
122  //BasicStructure &apply_sym(const SymOp &op); //Anna do this
123 
124  void generate_factor_group(SymGroup &factor_group, double map_tol) const;
125  void generate_factor_group_slow(SymGroup &factor_group, double map_tol) const;
126  void fg_converge(double small_tol, double large_tol, double increment);
127  void fg_converge(SymGroup &factor_group, double small_tol, double large_tol, double increment);
128 
129  SymGroupRepID generate_basis_permutation_representation(const MasterSymGroup &factor_group, bool verbose)const;
130 
131  void symmetrize(const SymGroup &relaxed_factors);
132  void symmetrize(const double &tolerance);
133 
134 
137  bool is_primitive(double prim_tol = TOL) const; //Donghee do this
138 
141  bool is_primitive(BasicStructure &new_prim, double prim_tol = TOL) const; //Donghee do this
142 
144  void fill_supercell(const BasicStructure &prim, double map_tol = TOL); //Ivy
145 
147  BasicStructure create_superstruc(const Lattice &scel_lat, double map_tol = TOL) const;
148 
149  //John G 230913
151  void generate_flowertrees_safe(const SiteOrbitree &in_tree, Array<SiteOrbitree> &out_trees);
152  void generate_flowertrees(const SiteOrbitree &in_tree, Array<SiteOrbitree> &out_trees);
153  //\John G 230913
154 
156  void map_superstruc_to_prim(BasicStructure &prim, const SymGroup &point_group);
157 
159  void merge_sites(double maxdist); //Only for same atom types
160 
161  //Array<Array<Array<double> > > get_NN_table(const double &maxr, GenericOrbitree<GenericCluster<CoordType> > &bouquet);
162  //Array<Array<Array<double> > > get_NN_table(const double &maxr);
163 
165  void add_vacuum_shift(BasicStructure &new_surface_struc, double vacuum_thickness, Eigen::Vector3d shift, COORD_TYPE mode) const;
166  void add_vacuum_shift(BasicStructure &new_surface_struc, double vacuum_thickness, Coordinate shift) const; //Because Anton thought a coordinate would be better
168  void add_vacuum(BasicStructure &new_surface_struc, double vacuum_thickness) const;
170  BasicStructure &operator+=(const Coordinate &shift);
171  BasicStructure &operator-=(const Coordinate &shift);
172 
175 
177  const Lattice &get_lattice() const {
178  return lattice;
179  };
180 
182  //CASM canonical input/output
183  virtual void read(std::istream &stream); //John do this
184 
186  void print_xyz(std::ostream &stream) const;
187  void print_cif(std::ostream &stream) const;
188 
189  jsonParser &to_json(jsonParser &json) const;
190 
191  // Assumes constructor CoordType::CoordType(Lattice) exists
192  void from_json(const jsonParser &json);
193  };
194 
195  template<typename CoordType>
196  BasicStructure<CoordType> operator*(const SymOp &LHS, const BasicStructure<CoordType> &RHS);
197 
198  template<typename CoordType>
199  BasicStructure<CoordType> operator*(const Lattice &LHS, const BasicStructure<CoordType> &RHS);
200 
201  //Translation operators -- not yet defined
202  template<typename CoordType>
203  BasicStructure<CoordType> operator+(const Coordinate &LHS, const BasicStructure<CoordType> &RHS);
204 
205  template<typename CoordType>
206  BasicStructure<CoordType> operator+(const BasicStructure<CoordType> &LHS, const Coordinate &RHS);
207 
208  template<typename CoordType>
209  jsonParser &to_json(const BasicStructure<CoordType> &basic, jsonParser &json);
210 
211 
212  // Assumes constructor CoordType::CoordType(Lattice) exists
213  template<typename CoordType>
214  void from_json(BasicStructure<CoordType> &basic, const jsonParser &json);
215 
218 };
219 
221 
222 #endif
void from_json(const jsonParser &json)
void print_cif(std::ostream &stream) const
void set_site_internals()
Associate each site with its basis index by setting its internal flags (asym_ind -> -1) ...
void map_superstruc_to_prim(BasicStructure &prim, const SymGroup &point_group)
Figures out which prim basis each superstructure basis corresponds to.
jsonParser & to_json(jsonParser &json) const
void fill_supercell(const BasicStructure &prim, double map_tol=TOL)
fill an empty structure with the basis of its corresponding primitive cell
void from_json(ClexDescription &desc, const jsonParser &json)
Type-safe ID object for communicating and accessing Symmetry representation info. ...
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
void add_vacuum(BasicStructure &new_surface_struc, double vacuum_thickness) const
Adds vacuum layer on top of ab plane.
bool is_primitive(double prim_tol=TOL) const
BasicStructure specifies the lattice and atomic basis of a crystal.
Definition: Cluster.hh:16
Unit Cell Coordinates.
BasicStructure create_superstruc(const Lattice &scel_lat, double map_tol=TOL) const
Shortcut routine to create a supercell structure and fill it with sites.
BasicStructure & operator+=(const Coordinate &shift)
Translates all atoms in cell.
void merge_sites(double maxdist)
If atoms are too close together, average their distance and make them one.
void generate_factor_group_slow(SymGroup &factor_group, double map_tol) const
Main CASM namespace.
Definition: complete.cpp:8
void generate_flowertrees(const SiteOrbitree &in_tree, Array< SiteOrbitree > &out_trees)
const Lattice & lattice() const
void main_print(std::ostream &stream, COORD_TYPE mode, bool version5, int option) const
const double TOL
void symmetrize(const SymGroup &relaxed_factors)
void set_lattice(const Lattice &lattice, COORD_TYPE mode)
const Lattice & get_lattice() const
Return a reference to the lattice.
void generate_factor_group(SymGroup &factor_group, double map_tol) const
apply a symmetry operation to the current structure (may change lattice vectors and order of basis at...
GenericCluster< CoordType > operator+(const GenericCluster< CoordType > &LHS, const Coordinate &RHS)
create translated cluster
Index find(const CoordType2 &test_site, double tol=TOL) const
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:33
BasicStructure(const Lattice &init_lat)
double tol
GenericOrbitree< SiteCluster > SiteOrbitree
Definition: Clex.hh:14
void print_xyz(std::ostream &stream) const
Output other formats.
BasisSet operator*(const SymOp &LHS, const BasisSet &RHS)
Definition: BasisSet.cc:1154
Represents cartesian and fractional coordinates.
Definition: Coordinate.hh:34
bool SD_flag
Specifies whether selectice dynamics is on or of for DFT calculations.
EigenIndex Index
For long integer indexing:
void add_vacuum_shift(BasicStructure &new_surface_struc, double vacuum_thickness, Eigen::Vector3d shift, COORD_TYPE mode) const
Add vacuum and shift c vector. The vacuum is always added parallel to c, and the shift vector should ...
void set_basis(Array< CoordType > basis_in)
Manually set the basis sites.
Array< CoordType > basis
Lattice vectors that specifies periodicity of the crystal.
Index max_possible_vacancies() const
Counts sites that allow vacancies.
void fg_converge(double small_tol, double large_tol, double increment)
std::string title
User-specified name of this Structure.
void generate_flowertrees_safe(const SiteOrbitree &in_tree, Array< SiteOrbitree > &out_trees)
Gets clusters of every size radiating from one site and saves them to a flowertree. A garland for each site is constructed.
BasicStructure & operator-=(const Coordinate &shift)
virtual void read(std::istream &stream)
Print intpolated images in seperate directries.
void copy_attributes_from(const BasicStructure &RHS)
CoordType get_site(const UnitCellCoord &ucc) const
UnitCellCoord get_unit_cell_coord(const CoordType2 &test_site, double tol=TOL) const
SymGroupRepID generate_basis_permutation_representation(const MasterSymGroup &factor_group, bool verbose) const
void within()
Translate all basis sites so that they are inside the unit cell.
virtual BasicStructure & operator=(const BasicStructure &RHS)