CASM  1.1.0
A Clusters Approach to Statistical Mechanics
PermuteIterator.cc
Go to the documentation of this file.
2 
3 #include <vector>
4 
9 #include "casm/external/Eigen/src/Core/Matrix.h"
11 
12 namespace CASM {
13 
15 
17  : m_sym_info(iter.m_sym_info),
18  m_trans_permute(&(m_sym_info->translation_permutations())),
19  m_factor_group_index(iter.m_factor_group_index),
20  m_translation_index(iter.m_translation_index) {}
21 
23  Index _factor_group_index,
24  Index _translation_index)
25  : m_sym_info(&_sym_info),
26  m_trans_permute(&(m_sym_info->translation_permutations())),
27  m_factor_group_index(_factor_group_index),
28  m_translation_index(_translation_index) {}
29 
31  swap(*this, iter);
32  return *this;
33 }
34 
36 const PermuteIterator &PermuteIterator::operator*() const { return *this; }
37 
39 const PermuteIterator *PermuteIterator::operator->() const { return this; }
40 
45 }
46 
49  return *m_sym_info;
50 }
51 
54  return sym_info().factor_group();
55 }
56 
62  return m_factor_group_index;
63 }
64 
70  return factor_group()[factor_group_index()].master_group_index();
71 }
72 
78 
81  return *(sym_info()
83  ->permutation());
84 }
85 
89 }
90 
94  return *sym_info().occ_symreps()[b][factor_group_index()];
95 }
96 
100  Index b) const {
101  return *sym_info().local_dof_symreps(_key)[b][factor_group_index()];
102 }
103 
107  DoFKey const &_key) const {
108  return *sym_info().global_dof_symrep(_key)[factor_group_index()];
109 }
110 
112  UnitCell translation_lattice_site =
114  SymOp lattice_translation_op = SymOp::translation(
115  make_superlattice_coordinate(translation_lattice_site,
116  this->sym_info().superlattice())
117  .cart());
118  return lattice_translation_op *
120 }
121 
124 }
125 
127  if (this->factor_group_index() == iter.factor_group_index()) {
128  return this->translation_index() < iter.translation_index();
129  }
130  return this->factor_group_index() < iter.factor_group_index();
131 }
132 
133 bool PermuteIterator::eq_impl(const PermuteIterator &iter) const {
134  if (m_sym_info == iter.m_sym_info &&
137  return true;
138  return false;
139 }
140 
141 // prefix ++PermuteIterator
144  if (m_translation_index == m_trans_permute->size()) {
147  }
148  return *this;
149 }
150 
151 // postfix PermuteIterator++
153  PermuteIterator cp(*this);
154  ++cp;
155  return cp;
156 }
157 
158 // prefix --PermuteIterator
160  if (m_translation_index == 0) {
163  }
165  return *this;
166 }
167 
168 // postfix PermuteIterator--
170  PermuteIterator cp(*this);
171  --cp;
172  return cp;
173 }
174 
176  PermuteIterator it(*this);
177  it.m_translation_index = 0;
179  return it;
180 }
181 
182 // Might be able to do this in a faster way, but would be much harder to
183 // understand
185  PermuteIterator it(*this);
186  // Finding the inverse factor_group operation is straightforward
188 
189  // Easiest way to get the new translation is just to compare the tau of the
190  // inverse of the 'total' sym_op (described by *this), to the inverse of the
191  // untranslated symop (described by
192  // m_fg_site_permutation_symrep.sym_op(it.m_factor_group_index)) Result is the
193  // portion of the inverse sym_op that needs to be described by a lattice point
194  // translation
195  Eigen::Vector3d translation_cart =
196  (sym_op().inverse().tau() -
197  factor_group()[it.factor_group_index()].tau());
198  UnitCell translation_uc = UnitCell::from_cartesian(
199  translation_cart, this->sym_info().prim_lattice());
201  this->sym_info().unitcell_index_converter()(translation_uc);
202 
203  return it;
204 }
205 
207  PermuteIterator it(*this);
208  // Finding the inverse factor_group operation is straightforward
211 
212  // Easiest way to get the new translation is just to compare the tau of the
213  // 'total' sym_op (described by (*this).sym_op()*RHS.sym_op()), to the
214  // untranslated symop product (described by
215  // m_fg_site_permutation_symrep.sym_op(it.factor_group_index())) Result is the
216  // portion of the product sym_op that needs to be described by a lattice point
217  // translation
218  Eigen::Vector3d translation_cart =
219  (sym_op() * RHS.sym_op()).tau() -
220  factor_group()[it.factor_group_index()].tau();
221  UnitCell translation_uc = UnitCell::from_cartesian(
222  translation_cart, this->sym_info().prim_lattice());
224  this->sym_info().unitcell_index_converter()(translation_uc);
225 
226  return it;
227 }
228 
235 template <typename PermuteIteratorIt>
236 SymGroup make_point_group(PermuteIteratorIt begin, PermuteIteratorIt end,
237  const Lattice &supercell_lattice) {
238  SymGroup result;
240  for (; begin != end; ++begin) {
241  Index f = begin->sym_op().index();
242  Index i;
243  for (i = 0; i < result.size(); ++i) {
244  if (f == result[i].index()) break;
245  }
246  if (i == result.size()) {
247  result.push_back((begin->sym_op()).no_trans());
248  // std::cout << "pushed back op " << result.back().index() << "\n";
249  }
250  }
251  result.sort();
252  return result;
253 }
254 
261 template <typename PermuteIteratorIt>
262 SymGroup make_sym_group(PermuteIteratorIt begin, PermuteIteratorIt end,
263  const Lattice &supercell_lattice) {
264  SymGroup sym_group;
265  sym_group.set_lattice(supercell_lattice);
266  while (begin != end) {
267  sym_group.push_back(begin->sym_op());
268  ++begin;
269  }
270  sym_group.sort();
271  return sym_group;
272 }
273 
281 template <typename PermuteIteratorIt>
282 std::unique_ptr<SymGroup> make_unique_sym_group(
283  PermuteIteratorIt begin, PermuteIteratorIt end,
284  const Lattice &supercell_lattice) {
285  auto sym_group_ptr = notstd::make_unique<SymGroup>();
286  sym_group_ptr->set_lattice(supercell_lattice);
287  while (begin != end) {
288  sym_group_ptr->push_back(begin->sym_op());
289  ++begin;
290  }
291  sym_group_ptr->sort();
292  return sym_group_ptr;
293 }
294 
296  const Lattice &supercell_lattice);
298  std::vector<PermuteIterator>::const_iterator begin,
299  std::vector<PermuteIterator>::const_iterator end,
300  const Lattice &supercell_lattice);
301 template SymGroup make_sym_group(std::vector<PermuteIterator>::iterator begin,
302  std::vector<PermuteIterator>::iterator end,
303  const Lattice &supercell_lattice);
304 
306  const Lattice &supercell_lattice);
308  std::vector<PermuteIterator>::const_iterator begin,
309  std::vector<PermuteIterator>::const_iterator end,
310  const Lattice &supercell_lattice);
311 template SymGroup make_point_group(std::vector<PermuteIterator>::iterator begin,
312  std::vector<PermuteIterator>::iterator end,
313  const Lattice &supercell_lattice);
314 
320 }
321 
324  std::set<Index> const &site_indices) {
325  // Applying the permutation indicated by `permute_it` moves the value from
326  // site index `permute_it.permute_ind(s)` to site index `s`, for each `s` in
327  // the set. Therefore, if none of `permute_it.permute_ind(s)` are outside the
328  // set `site_indices` the sites are invariant.
329 
330  return std::none_of(site_indices.begin(), site_indices.end(), [&](Index s) {
331  return site_indices.count(permute_it.permute_ind(s)) == 0;
332  });
333 }
334 
336  json.put_obj();
337  json["factgrp"] = it.factor_group_index();
338  json["trans"] = it.translation_index();
339  return json;
340 }
341 
343  const jsonParser &json, const SupercellSymInfo &scel_info) {
344  return scel_info.permute_it(json["factgrp"].get<Index>(),
345  json["trans"].get<Index>());
346 }
347 
348 namespace adapter {
349 
351  PermuteIterator const &adaptable) const {
352  return adaptable.sym_op();
353 }
354 } // namespace adapter
355 
356 } // namespace CASM
std::set< std::string > & s
Index prim_factor_group_index() const
Return the prim factor group index.
bool operator<(const PermuteIterator &iter) const
SupercellSymInfo const * m_sym_info
SupercellSymInfo const & sym_info() const
Reference the SupercellSymInfo containing the operations being pointed at.
PermuteIterator begin_next_fg_op() const
SymOpRepresentation const & global_dof_rep(DoFKey const &_key) const
Index factor_group_index() const
Return the supercell factor group index.
const Permutation & translation_permute() const
Return the translation permutation being pointed at.
SymOpRepresentation const & local_dof_rep(DoFKey const &_key, Index b) const
bool eq_impl(const PermuteIterator &iter) const
PermuteIterator inverse() const
PermuteIterator & operator--()
Index translation_index() const
Return the index into the supercell translation permutations.
const PermuteIterator * operator->() const
Returns a pointer to this.
const Permutation & factor_group_permute() const
Return the factor group permutation being pointed at.
SymOpRepresentation const & occ_rep(Index b) const
const PermuteIterator & operator*() const
Returns a reference to this.
friend void swap(PermuteIterator &a, PermuteIterator &b)
Permutation combined_permute() const
Index permute_ind(Index i) const
PermuteIterator & operator=(PermuteIterator iter)
PermuteIterator & operator++()
std::vector< Permutation > const * m_trans_permute
SymGroup const & factor_group() const
Reference the SupercellSymInfo::factor_group()
A class that collects all symmetry information for for performing symmetry transformations on the sit...
SymGroup const & factor_group() const
Subgroup of primitive-cell factor group operations that leave supercell lattice invariant.
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
permute_const_iterator permute_it(Index supercell_factor_group_index, Index translation_index) const
std::map< DoFKey, SublatSymReps > const & local_dof_symreps() const
Const reference local DoF matrix representations of the supercell's factor group.
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
Index ind_inverse(Index i) const
Get index of operation that is inverse of operation at(i)
Definition: SymGroup.cc:1540
virtual void sort()
Sort SymOp in the SymGroup.
Definition: SymGroup.cc:2016
Index ind_prod(Index i, Index j) const
Get index of operation that is result of multiplication of at(i)*at(j)
Definition: SymGroup.cc:1550
virtual void push_back(const SymOp &new_op)
Definition: SymGroup.cc:818
void set_lattice(const Lattice &new_lat)
Lattice used for periodic comparisons (for instance, to generate multiplcation table)
Definition: SymGroup.cc:814
Permutation const * permutation(Index i) const
pointer to Permutation corresponding to SymOpRepresentation at entry 'i' of this SymGroupRep Returns ...
Definition: SymGroupRep.cc:157
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
Definition: SymOp.hh:28
SymOp inverse() const
get the inverse of this SymOp
Definition: SymOp.cc:87
const vector_type & tau() const
Const access of the cartesian translation vector, 'tau'.
Definition: SymOp.hh:63
static SymOp translation(const Eigen::Ref< const vector_type > &_tau)
static method to create operation that describes pure translation
Definition: SymOp.hh:34
SymOpRepresentation is the base class for anything describes a symmetry operation.
jsonParser & put_obj()
Puts new empty JSON object.
Definition: jsonParser.hh:354
Unit Cell Indices.
static UnitCell from_cartesian(const Eigen::Vector3d &cartesian_coord, const Lattice &tiling_unit)
VectorXdSupercellSymInfoFormatter supercell_lattice()
SupercellSymInfoFormatter< jsonParser > translation_permutations()
Coordinate make_superlattice_coordinate(const UnitCell &ijk, const Superlattice &superlattice)
Main CASM namespace.
Definition: APICommand.hh:8
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
SymGroup make_sym_group(const PermuteIteratorContainer &container, const Lattice &supercell_lattice)
Returns a SymGroup generated from a container of PermuteIterator.
void swap(ConfigDoF &A, ConfigDoF &B)
Definition: ConfigDoF.cc:260
bool site_indices_are_invariant(PermuteIterator const &permute_it, std::set< Index > const &site_indices)
Return true if the permutation does not mix given sites and other sites.
std::string DoFKey
Definition: DoFDecl.hh:7
std::unique_ptr< SymGroup > make_unique_sym_group(const PermuteIteratorContainer &container, const Lattice &supercell_lattice)
Returns a std::unique_ptr<SymGroup> generated from a container of PermuteIterator.
void swap(PermuteIterator &a, PermuteIterator &b)
INDEX_TYPE Index
For long integer indexing:
Definition: definitions.hh:39
SymGroup make_point_group(const PermuteIteratorContainer &container, const Lattice &supercell_lattice)
Returns a SymGroup generated from a container of PermuteIterator.
static ReturnType from_json(const jsonParser &json)
Default from_json is equivalent to.
Definition: jsonParser.hh:551