CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
SymGroupRep.hh
Go to the documentation of this file.
1 #ifndef SYMGROUPREP_HH
2 #define SYMGROUPREP_HH
3 
4 #include <iostream>
5 #include <string>
6 #include <iomanip>
7 
10 #include "casm/symmetry/SymOp.hh"
12 
13 
14 namespace CASM {
15  class SymGroupRep;
16  class MasterSymGroup;
17  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
18  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
19 
20  class SymGroupRepHandle;
21 
26  class SymGroupRep : public Array<SymOpRepresentation *> {
31  public:
34 
37 
39  Array<SymOpRepresentation * >(_size, nullptr),
40  m_master_group(nullptr) { }
41 
42  SymGroupRep(const SymGroup &_head, SymGroupRepID _rep_ID = SymGroupRepID()) :
43  m_rep_ID(_rep_ID), m_master_group(&_head.master_group()) {
44  if(has_valid_master())
45  resize(master_group().size(), nullptr);
46  else {
47  std::cerr << "WARNING: SymGroupRep initialized without valid master_group! Bad things may happen!\n";
48  assert(0);
49  }
50  }
51 
52  SymGroupRep(const SymGroupRep &RHS);
53 
54  ~SymGroupRep();
55 
56  //void push_back(const SymOpRepresentation &new_op);
57 
62  void set_rep(const SymOpRepresentation &base_rep, const SymOpRepresentation &new_rep) const;
63  void set_rep(const SymOp &base_rep, const SymOpRepresentation &new_rep) const;
64  void set_rep(Index op_index, const SymOpRepresentation &new_rep) const;
65 
66  const MasterSymGroup &master_group() const {
67  assert(m_master_group);
68  return *m_master_group;
69  }
70 
71  bool has_valid_master() const {
72  return m_master_group != nullptr;
73  }
74  void set_master_group(const MasterSymGroup &master, const SymGroupRepID &_rep_ID);
75 
78 
79  SymGroupRep &operator=(const SymGroupRep &RHS);
80 
81  void clear();
82 
83  SymGroupRep *copy()const {
84  return new SymGroupRep(*this);
85  }
86 
87  //John G 050513
88 
89  void print_permutation(std::ostream &stream) const;
90  void print_MatrixXd(std::ostream &stream) const;
91  void print_MatrixXd(std::ostream &stream, const SymGroup &subgroup) const;
92 
93  //\John G
94 
95  // block_shape_matrix is sum of squares of each (i,j) matrix element over all operations in SymGroupRep
96  // It reveals the block_diagonalization of the symgrouprep
98  Eigen::MatrixXd block_shape_matrix(const SymGroup &subgroup) const;
99 
100  // counts number of blocks in the block_shape_matrix
101  Index num_blocks() const;
102  Index num_blocks(const SymGroup &subgroup) const;
103 
104  Eigen::MatrixXd const *get_MatrixXd(Index i) const;
105  Eigen::MatrixXd const *get_MatrixXd(const SymOpRepresentation &) const;
106  Permutation const *get_permutation(Index i) const;
107  Permutation const *get_permutation(const SymOpRepresentation &) const;
108 
109  // adds a copy of '_pushed' temporary workaround for 'masterless' SymGroupReps; to be deprecated
110  // will throw if valid master exists.
111  void push_back_copy(const SymOpRepresentation &_pushed);
112 
114  return m_rep_ID;
115  }
116 
117  std::vector<Eigen::MatrixXd> irreducible_wedges(const SymGroup &head_group, std::vector<Index> &multiplicities)const;
120 
122  ReturnArray<Index> num_each_irrep(const SymGroup &sub_group, bool verbose = false) const;
123  ReturnArray<Index> num_each_real_irrep(const SymGroup &subgroup, bool verbose = false) const;
124 
125  ReturnArray<SymGroupRepID> get_irrep_IDs(const SymGroup &subgroup) const;
126 
127  bool is_irrep() const;
128  bool is_irrep(const SymGroup &head_group) const;
129 
133  SymGroupRep coord_transformed_copy(const Eigen::MatrixXd &trans_mat) const;
134 
137  SymGroupRep symmetry_adapted_copy(const SymGroup &head_group) const {
138  return coord_transformed_copy(get_irrep_trans_mat(head_group));
139  }
140 
144  Eigen::MatrixXd get_irrep_trans_mat(const SymGroup &head_group) const;
149  Eigen::MatrixXd get_irrep_trans_mat(const SymGroup &head_group, std::vector<std::vector<Index> > &subspaces) const;
150 
151  Eigen::MatrixXd get_irrep_trans_mat_blind(const SymGroup &head_group) const;
153 
154  jsonParser &to_json(jsonParser &json) const;
155 
156  // If 'm_master_group' is not nullptr, should be initialized accordingly
157  void from_json(const jsonParser &json);
158  private:
159 
162 
165  // don't use default constructor
167  std::cerr << "Cannot perform default construction of SymGroupRep.\nExiting...\n";
168  exit(1);
169  }
170 
173 
175  SymGroupRep(const MasterSymGroup *_home) : m_master_group(_home) { }
176 
177  void calc_new_irreps(int max_iter = 1000) const;
178  void calc_new_irreps(const SymGroup &sub_group, int max_iter = 1000) const;
179 
181 
187 
188  };
189 
190  SymGroupRep subset_permutation_rep(const SymGroupRep &permute_rep, const Array<Index>::X2 &subsets);
191  SymGroupRep permuted_direct_sum_rep(const SymGroupRep &permute_rep, const Array<SymGroupRep const *> &sum_reps);
192  SymGroupRep kron_rep(const SymGroupRep &LHS, const SymGroupRep &RHS);
193 
194  jsonParser &to_json(const SymGroupRep &rep, jsonParser &json);
195 
196  // If 'm_home_group' is not nullptr, should be initialized accordingly
197  void from_json(SymGroupRep &rep, const jsonParser &json);
198 
199  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
200  //~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
201 
204  class SymGroupRepHandle { // <-- typedefed as SymGroupRep::RemoteHandle
207  public:
209  m_group_rep(nullptr) {}
210 
211  SymGroupRepHandle(const SymGroup &head_group, SymGroupRepID symrep_ID) :
212  m_group_rep(nullptr), m_subgroup_op_inds(head_group.op_indices()) {
213  if(head_group.size() == 0 || !head_group[0].has_valid_master()) {
214  std::cerr << "CRITICAL ERROR: Requested representation of an improperly initialized SymGroup.\n"
215  << " Exiting...\n";
216  assert(0);
217  }
218  assert(!symrep_ID.empty());
219  m_group_rep = &(head_group[0].master_group().representation(symrep_ID));
220  assert(m_group_rep);
221  }
222 
224  Index size() const {
225  return m_subgroup_op_inds.size();
226  }
227 
229  Index dim() const {
230  return (m_group_rep->get_MatrixXd(m_subgroup_op_inds[0]))->cols();
231  }
232 
233  SymGroupRep const *rep_ptr() const {
234  return m_group_rep;
235  }
236 
237  SymGroupRep const *operator->() const {
238  return rep_ptr();
239  }
240 
242  return m_group_rep->at(m_subgroup_op_inds[i]);
243  }
244 
245  const SymOp &sym_op(Index i) const {
246  return (m_group_rep->master_group())[m_subgroup_op_inds[i]];
247  }
248 
250  return m_subgroup_op_inds.find((m_group_rep->master_group()).ind_inverse(m_subgroup_op_inds[i]));
251  }
252 
253  Index ind_prod(Index i, Index j) const {
254  return m_subgroup_op_inds.find((m_group_rep->master_group()).ind_prod(m_subgroup_op_inds[i], m_subgroup_op_inds[j]));
255  }
256 
257  bool operator==(const SymGroupRepHandle &RHS) const {
258  return m_group_rep == RHS.m_group_rep && m_subgroup_op_inds == RHS.m_subgroup_op_inds;
259  }
260 
261  void set_rep(const SymGroup &head_group, SymGroupRepID symrep_ID) {
262  m_subgroup_op_inds = head_group.op_indices();
263  assert(head_group.size() && !symrep_ID.empty() && head_group[0].has_valid_master());
264  m_group_rep = &(head_group[0].master_group().representation(symrep_ID));
265  assert(m_group_rep);
266 
267  }
268  };
269 
271 }
272 
273 #endif
Eigen::MatrixXd MatrixXd
jsonParser & to_json(jsonParser &json) const
ReturnArray< Eigen::MatrixXd > get_projection_operators() const
Index size() const
Size of the associate SymGroup.
Definition: SymGroupRep.hh:224
std::vector< Eigen::MatrixXd > irreducible_wedges(const SymGroup &head_group, std::vector< Index > &multiplicities) const
Definition: SymGroupRep.cc:442
Index ind_prod(Index i, Index j) const
Definition: SymGroupRep.hh:253
void from_json(ClexDescription &desc, const jsonParser &json)
const MasterSymGroup & master_group() const
Definition: SymGroup.hh:56
Type-safe ID object for communicating and accessing Symmetry representation info. ...
bool empty() const
Returns true if SymGroupRepID has not been initialized with valid group_index or rep_index.
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
ReturnArray< Index > num_each_irrep() const
Definition: SymGroupRep.cc:695
SymGroupRep subset_permutation_rep(const SymGroupRep &permute_rep, const Array< Index >::X2 &subsets)
ReturnArray< Array< Eigen::VectorXd > > _calc_special_irrep_directions(const SymGroup &subgroup) const
Definition: SymGroupRep.cc:471
SymGroupRep symmetry_adapted_copy(const SymGroup &head_group) const
Definition: SymGroupRep.hh:137
SymGroupRep const * operator->() const
Definition: SymGroupRep.hh:237
Eigen::MatrixXd block_shape_matrix() const
Definition: SymGroupRep.cc:180
Main CASM namespace.
Definition: complete.cpp:8
void from_json(const jsonParser &json)
ReturnArray< SymGroupRepID > get_irrep_IDs(const SymGroup &subgroup) const
Definition: SymGroupRep.cc:747
bool is_irrep() const
Definition: SymGroupRep.cc:981
SymGroupRep(const SymGroup &_head, SymGroupRepID _rep_ID=SymGroupRepID())
Definition: SymGroupRep.hh:42
Eigen::MatrixXd get_irrep_trans_mat(const SymGroup &head_group) const
const SymOp & sym_op(Index i) const
Definition: SymGroupRep.hh:245
void print_MatrixXd(std::ostream &stream) const
Definition: SymGroupRep.cc:152
SymGroupRepHandle(const SymGroup &head_group, SymGroupRepID symrep_ID)
Definition: SymGroupRep.hh:211
SymGroupRep kron_rep(const SymGroupRep &LHS, const SymGroupRep &RHS)
SymGroup is a collection of symmetry operations that satisfy the group property The symmetry operatio...
Definition: SymGroup.hh:33
Array< X1 > X2
Definition: Array.hh:64
void calc_new_irreps(int max_iter=1000) const
Definition: SymGroupRep.cc:772
MasterSymGroup const * m_master_group
Pointer to the home_group that generated this SymGroupRep.
Definition: SymGroupRep.hh:164
bool has_valid_master() const
Definition: SymGroupRep.hh:71
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
Definition: SymOp.hh:28
ReturnArray< Array< Eigen::MatrixXd > > calc_special_subspaces(const SymGroup &subgroup) const
Definition: SymGroupRep.cc:568
typename multivector_impl::multivector_tmp< T, N >::type X
Definition: multivector.hh:28
Eigen::MatrixXd const * get_MatrixXd(Index i) const
Definition: SymGroupRep.cc:253
SymGroupRep(const MasterSymGroup *_home)
Pointer version of constructor is private for internal construction of master-less representations...
Definition: SymGroupRep.hh:175
SymGroupRep const * m_group_rep
Definition: SymGroupRep.hh:205
EigenIndex Index
For long integer indexing:
Index ind_inverse(Index i) const
Definition: SymGroupRep.hh:249
void print_permutation(std::ostream &stream) const
Definition: SymGroupRep.cc:138
Index num_blocks() const
Definition: SymGroupRep.cc:207
bool operator==(const SymGroupRepHandle &RHS) const
Definition: SymGroupRep.hh:257
SymGroupRep const * rep_ptr() const
Definition: SymGroupRep.hh:233
Index dim() const
Matrix dimension of representation.
Definition: SymGroupRep.hh:229
SymOpRepresentation is the base class for anything describes a symmetry operation.
SymGroupRepID add_copy_to_master() const
Adds copy of this representation its home_group.
Definition: SymGroupRep.cc:123
SymGroupRep & operator=(const SymGroupRep &RHS)
Definition: SymGroupRep.cc:27
Index find(const T &test_elem) const
Definition: Array.hh:707
void push_back_copy(const SymOpRepresentation &_pushed)
Definition: SymGroupRep.cc:244
ReturnArray< Index > num_each_real_irrep(const SymGroup &subgroup, bool verbose=false) const
Definition: SymGroupRep.cc:668
const MasterSymGroup & master_group() const
Definition: SymGroupRep.hh:66
Permutation const * get_permutation(Index i) const
Definition: SymGroupRep.cc:265
void set_rep(const SymGroup &head_group, SymGroupRepID symrep_ID)
Definition: SymGroupRep.hh:261
SymGroupRep permuted_direct_sum_rep(const SymGroupRep &permute_rep, const Array< SymGroupRep const * > &sum_reps)
SymGroupRep coord_transformed_copy(const Eigen::MatrixXd &trans_mat) const
Definition: SymGroupRep.cc:280
SymGroupRepHandle RemoteHandle
Definition: SymGroupRep.hh:32
Eigen::MatrixXd _symmetrized_irrep_trans_mat(const SymGroup &subgroup) const
Definition: SymGroupRep.cc:304
SymGroupRep is an alternative representation of a SymGroup for something other than real space...
Definition: SymGroupRep.hh:30
SymGroupRepID get_ID() const
Definition: SymGroupRep.hh:113
multivector< Eigen::VectorXd >::X< 3 > calc_special_total_directions(const SymGroup &subgroup) const
Definition: SymGroupRep.cc:415
Eigen::MatrixXd get_irrep_trans_mat_blind(const SymGroup &head_group) const
void set_master_group(const MasterSymGroup &master, const SymGroupRepID &_rep_ID)
Definition: SymGroupRep.cc:43
SymGroupRepID m_rep_ID
rep_ID is unique identifier of a specific SymGroupRep instantiation
Definition: SymGroupRep.hh:161
SymGroupRep(SymGroupRep::NullInitializer init, Index _size)
Definition: SymGroupRep.hh:38
SymOpRepresentation const * operator[](Index i) const
Definition: SymGroupRep.hh:241
Array< Index > m_subgroup_op_inds
Definition: SymGroupRep.hh:206
SymGroupRep * copy() const
Definition: SymGroupRep.hh:83
Basic std::vector like container (deprecated)
ReturnArray< Index > op_indices() const
Definition: SymGroup.cc:474
SymGroupRep const & representation(SymGroupRepID i) const
Const access of alternate Representations of a SymGroup.
Definition: SymGroup.cc:375
void set_rep(const SymOpRepresentation &base_rep, const SymOpRepresentation &new_rep) const
Definition: SymGroupRep.cc:82