CASM  1.1.0
A Clusters Approach to Statistical Mechanics
SymTools.hh
Go to the documentation of this file.
1 #ifndef XTALSYMTOOLS_HH
2 #define XTALSYMTOOLS_HH
3 
4 #include <vector>
5 
10 
11 namespace CASM {
12 namespace sym {
13 // TODO: These could be template specializations of what's in
14 // symmetry/SymTools.hh but I'm not sure how we're generalizing the
15 // apply/copy_apply stuff yet
16 
18 xtal::Lattice &apply(const xtal::SymOp &op, xtal::Lattice &lat);
19 
21 xtal::Lattice copy_apply(const xtal::SymOp &op, xtal::Lattice lat_copy);
22 
23 template <typename ExternSymOp>
24 xtal::Lattice copy_apply(const ExternSymOp &op, const xtal::Lattice &lat) {
26 }
27 } // namespace sym
28 
29 namespace xtal {
30 
32 std::vector<Index> invariant_subgroup_indices(const Lattice &lat,
33  SymOpVector const &super_grp);
34 
35 template <typename ExternSymOpVector>
36 std::vector<Index> invariant_subgroup_indices(
37  const Lattice &lat, ExternSymOpVector const &super_grp) {
40 }
41 
44 template <typename OutputIt>
45 OutputIt invariant_subgroup_indices(const Lattice &lat,
46  const std::vector<SymOp> &super_group,
47  OutputIt result) {
48  IsPointGroupOp is_equiv(lat);
49 
50  Index ix = 0;
51  for (auto it = super_group.begin(); it != super_group.end(); ++it) {
52  if (is_equiv(*it)) {
53  *result = ix;
54  ++result;
55  }
56  ++ix;
57  }
58  return result;
59 }
60 
61 // TODO
62 // Unimplemented. Is this the same as cart2frac and frac2cart?
64 /* Eigen::Matrix3i symmetry_matrix_to_frac(const Lattice &lat, const
65  * Eigen::Matrix3d &cart_matrix); */
67 /* Eigen::Matrix3d symmetry_matrix_to_cart(const Lattice &lat, const
68  * Eigen::Matrix3i &cart_matrix); */
69 
72 Lattice symmetrize(const Lattice &lat,
73  const std::vector<SymOp> &enforced_group);
74 
75 // TODO
76 // Why does this routine take a tolerance, if the lattice itself has a
77 // tolerance? Should we get rid of the tolerance inside of lattice?
81 Lattice symmetrize(const Lattice &lat, double point_group_tolerance);
82 
86 /* Lattice symmetrized_with_fractional(const Lattice& lat, const
87  * std::vector<Eigen::Matrix3i> &fractional_point_group);
88  */
89 
94 std::vector<SymOp> make_point_group(Lattice const &_lat);
95 
100 std::vector<SymOp> make_point_group(Lattice const &_lat, double _tol);
101 
102 //************************************************************************************************************************//
103 
109 template <typename Object, typename OpIterator>
110 std::pair<OpIterator, Eigen::Matrix3d> is_equivalent_superlattice(
111  const Object &scel, const Object &unit, OpIterator begin, OpIterator end,
112  double tol) {
113  std::pair<bool, Eigen::Matrix3d> res;
114  for (auto it = begin; it != end; ++it) {
115  res = is_superlattice(scel, sym::copy_apply(*it, unit), tol);
116  if (res.first) {
117  return std::make_pair(it, res.second);
118  }
119  }
120  return std::make_pair(end, res.second);
121 }
122 
129 template <typename LatIterator, typename SymOpIterator>
130 Lattice make_equivalent_superduperlattice(LatIterator begin, LatIterator end,
131  SymOpIterator op_begin,
132  SymOpIterator op_end) {
133  Lattice best = *begin;
134  for (auto it = ++begin; it != end; ++it) {
135  Lattice tmp_best = make_superduperlattice(best, *it);
136  for (auto op_it = op_begin; op_it != op_end; ++op_it) {
137  Lattice test = make_superduperlattice(best, sym::copy_apply(*op_it, *it));
138  if (std::abs(volume(test)) < std::abs(volume(tmp_best))) {
139  tmp_best = test;
140  }
141  }
142  best = tmp_best;
143  }
144  return best;
145 }
146 
147 } // namespace xtal
148 
149 } // namespace CASM
150 
151 namespace CASM {
152 namespace xtal {
153 
154 template <typename StructureType, typename ExternSymOpVector>
155 StructureType symmetrize(const StructureType &structure,
156  const ExternSymOpVector &enforced_group) {
157  return xtal::symmetrize(
158  structure,
160 }
161 } // namespace xtal
162 } // namespace CASM
163 
164 #endif
Checks if operations are point group operations.
double volume(const Lattice &lat)
Returns the volume of a Lattice.
Definition: Lattice.hh:281
std::pair< bool, Eigen::Matrix3d > is_superlattice(const Lattice &scel, const Lattice &unit, double tol)
Check if scel is a superlattice of unitcell unit and some integer transformation matrix T.
Definition: Lattice.cc:836
Lattice make_superduperlattice(const Lattice &lat1, const Lattice &lat2)
returns Lattice that is smallest possible superlattice of both input Lattice
Definition: Lattice.cc:783
ConfigIO::GenericConfigFormatter< jsonParser > structure()
Definition: ConfigIO.cc:766
xtal::Coordinate copy_apply(const xtal::SymOp &op, xtal::Coordinate coord)
Copy and apply SymOp to a Coordinate.
Definition: Coordinate.cc:354
xtal::Coordinate & apply(const xtal::SymOp &op, xtal::Coordinate &coord)
apply SymOp to a Coordinate
Definition: Coordinate.cc:347
BasicStructure symmetrize(const BasicStructure &structure, const std::vector< SymOp > &enforced_group)
std::vector< Index > invariant_subgroup_indices(const Lattice &lat, SymOpVector const &super_grp)
Construct indices of the subgroup that leaves a lattice unchanged.
Lattice make_equivalent_superduperlattice(LatIterator begin, LatIterator end, SymOpIterator op_begin, SymOpIterator op_end)
returns Lattice that is smallest possible superlattice of all input Lattice
Definition: SymTools.hh:130
std::pair< OpIterator, Eigen::Matrix3d > is_equivalent_superlattice(const Object &scel, const Object &unit, OpIterator begin, OpIterator end, double tol)
Definition: SymTools.hh:110
std::vector< SymOp > make_point_group(Lattice const &_lat)
Populate.
Definition: SymTools.cc:113
std::vector< SymOp > SymOpVector
Main CASM namespace.
Definition: APICommand.hh:8
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39