CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
SymOp.cc
Go to the documentation of this file.
1 #include "casm/symmetry/SymOp.hh"
2 
7 namespace CASM {
8 
9  const double &SymOp::map_error() const {
10  return m_map_error;
11  }
12 
13  //*******************************************************************************************
14 
15  void SymOp::set_map_error(const double &value) {
16  m_map_error = value;
17  return;
18  }
19 
20 
21  //*******************************************************************************************
22 
23  void SymOp::set_index(const MasterSymGroup &new_group, Index new_index) {
24  if((valid_index(new_index) && new_index < new_group.size())
25  && (this == &(new_group[new_index]) ||
26  almost_equal(matrix(), new_group[new_index].matrix()))) {
27  m_master_group = &new_group;
28  m_op_index = new_index;
29  }
30  else {
31  m_master_group = &new_group;
32  //m_master_group = NULL;
33  m_op_index = -1;
34  }
35  }
36 
37  //*******************************************************************************************
38 
39  SymOp SymOp::operator*(const SymOp &RHS) const {
40  SymOp t_op(matrix() * RHS.matrix(),
41  tau() + matrix() * RHS.tau());
42 
45  }
46  else if(RHS.m_master_group && !m_master_group && is_identity()) {
47  t_op.set_index(*RHS.m_master_group, RHS.index());
48  }
49  else if(m_master_group && !RHS.m_master_group && RHS.is_identity()) {
50  t_op.set_index(master_group(), index());
51  }
52  //The following blocks caused problems at some point (mainly for non-primitive structures)
53  //else if(is_identity() && RHS.is_identity()) {
54  //t_op.m_op_index = 0;
55  //}
56  //else{
57  //std::cout << "This symmetry is " << symmetry << " with head " << m_master_group << " and RHS symmetry is " << RHS.symmetry << " with head " << RHS.m_master_group << "\n";
58  //}
59 
60  return t_op;
61 
62  }
63 
64  //*******************************************************************************************
65 
66  SymOp &SymOp::operator+=(const Eigen::Ref<const SymOp::vector_type> &RHS) {
67  m_tau += RHS - matrix() * RHS;
68  return (*this);
69  }
70 
71  //*******************************************************************************************
72 
73  SymOp &SymOp::operator-=(const Eigen::Ref<const SymOp::vector_type> &RHS) {
74  m_tau -= RHS - matrix() * RHS;
75  return (*this);
76  }
77 
78  //*******************************************************************************************
79  // SymOp matrix is unitary, so inverse is equivalent to transpose.
80  // To do inverse of translation, you must perform
81  // inverse matrix operaton on translation and subtract
83  SymOp t_op(matrix().transpose(),
84  -(matrix().transpose() * tau()));
85  if(m_master_group) {
87  }
88  else if(is_identity()) {
89  t_op.m_op_index = 0;
90  }
91 
92  return t_op;
93  }
94 
95  //*******************************************************************************************
96 
98  return SymOp(matrix(), vector_type::Zero(), map_error(), index(), m_master_group);
99  }
100 
101  //*******************************************************************************************
102 
103  bool SymOp::operator==(const SymOp &RHS) const {
104  return
105  almost_equal(matrix(), RHS.matrix()) &&
106  almost_equal(tau(), RHS.tau());
107  };
108 
109  //*******************************************************************************************
110 
112  (*this) = op * (*this) * (op.inverse());
113  return *this;
114  }
115 
116  //*******************************************************************************************
117 
118  void SymOp::print(std::ostream &stream, const Eigen::Ref<const SymOp::matrix_type> &c2f_mat) const {
119 
120  int tprec = stream.precision();
121  std::ios::fmtflags tflags = stream.flags();
122 
123  stream.precision(3);
124 
125  stream.flags(std::ios::left);
126  stream << std::setw(53) << "Symmetry Operation Matrix" << "Shift \n"; //SOM has 25 chars, width of those 3 columns are 14 each, so 42 width. Shift width is 22, so spacing of 9-11 extra characters, so add 5 more to get in the middle
127 
128  stream.flags(std::ios::showpoint | std::ios::fixed | std::ios::right);
129  stream.precision(9);
130  matrix_type tmat(c2f_mat * matrix()*c2f_mat.inverse());
131  vector_type ttau(c2f_mat * tau());
132  for(int i = 0; i < 3; i++) {
133  //Print each row of the symmetry matrix separately
134  for(int j = 0; j < 3; j++) {
135  stream << std::setw(14) << tmat(i, j);
136  }
137 
138  stream << std::setw(22) << ttau(i) << "\n";
139  }
140 
141  stream.precision(tprec);
142  stream.flags(tflags);
143  return;
144  }
145 
146  //*******************************************************************************************
147 
149  json.put_obj();
150 
151  // Members not included:
152  //
153  // From SymOpRepresentation:
154  // MasterSymGroup const master_group();
155  //
156  // From SymOp:
157  // Lattice const *home;
158  // Array<SymOpRepresentation *> representation;
159 
160  json["SymOpRep_type"] = "SymOp";
161 
162  json["op_index"] = m_op_index;
163  json["rep_ID"] = m_rep_ID;
164 
165  // mutable SymOp::matrix_type symmetry_mat[2];
166  json["symmetry_mat"] = matrix();
167 
168  // mutable Coordinate tau_vec;
169  to_json_array(tau(), json["tau"]);
170 
171  // double map_error;
172  json["map_error"] = map_error();
173 
174  return json;
175  }
176 
177  //*******************************************************************************************
178 
179  void SymOp::from_json(const jsonParser &json) {
180  try {
181  //std::cout<<"Inside of SymOp::from_json"<<std::endl;
182  //std::cout<<"Reading in op_index"<<std::endl;
183  CASM::from_json(m_op_index, json["op_index"]);
184  //std::cout<<"Reading in rep_id"<<std::endl;
185  CASM::from_json(m_rep_ID, json["rep_ID"]);
186 
187  // mutable SymOp::matrix_type symmetry_mat[2];
188  //std::cout<<"Reading in symmetry_mat"<<std::endl;
189  CASM::from_json(m_mat, json["symmetry_mat"]);
190 
191  // mutable Coordinate tau_vec;
192  //std::cout<<"Reading in tau_vec"<<std::endl;
193  if(json.contains("tau"))
194  CASM::from_json(m_tau, json["tau"]);
195 
196  //std::cout<<"Reading in map_error"<<std::endl;
197  // double map_error;
198  CASM::from_json(m_map_error, json["map_error"]);
199  //std::cout<<"Done Reading in the SymOp"<<std::endl;
200  }
201  catch(...) {
203  throw;
204  }
205  }
206 
207  //*******************************************************************************************
208 
209  jsonParser &to_json(const SymOp &sym, jsonParser &json) {
210  return sym.to_json(json);
211  }
212 
213  //*******************************************************************************************
214  void from_json(SymOp &sym, const jsonParser &json) {
215  try {
216  sym.from_json(json);
217  }
218  catch(...) {
220  throw;
221  }
222  }
223 
224 
225 }
SymOp(const Eigen::Ref< const matrix_type > &_mat=matrix_type::Identity(), const Eigen::Ref< const vector_type > &_tau=vector_type::Zero(), double _map_error=TOL)
Definition: SymOp.hh:41
SymOp operator*(const SymOp &RHS) const
get a new SymOp that is equivalent to subsequent application of both SymOps
Definition: SymOp.cc:39
void from_json(ClexDescription &desc, const jsonParser &json)
Index size() const
Definition: Array.hh:145
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
void set_index(const MasterSymGroup &new_group, Index new_index)
Definition: SymOp.cc:23
bool operator==(const SymOp &RHS) const
Check equality of SymOps, (matrix and translation). Does not necessarily return true for translationa...
Definition: SymOp.cc:103
MasterSymGroup const * m_master_group
Pointer to the MasterSymGroup where prototype of this SymOp lives.
Index index() const
Index of this operation within the master_group.
SymOp & operator-=(const Eigen::Ref< const vector_type > &RHS)
Definition: SymOp.cc:73
Main CASM namespace.
Definition: complete.cpp:8
SymGroupRepID m_rep_ID
ID of this representation within the master_group. Default is uninitialized.
Index m_op_index
Index into MasterSymGroup that specifies the operation.
bool is_identity() const
returns true if matrix part of operation is identity
Definition: SymOp.hh:69
void set_map_error(const double &value)
Definition: SymOp.cc:15
jsonParser & to_json(jsonParser &json) const override
Definition: SymOp.cc:148
matrix_type m_mat
Definition: SymOp.hh:158
const matrix_type & matrix() const
Const access of entire cartesian symmetry matrix.
Definition: SymOp.hh:57
Eigen::Vector3d vector_type
Definition: SymOp.hh:32
void from_json(const jsonParser &json) override
Definition: SymOp.cc:179
const double & map_error() const
Allows access to the map_error.
Definition: SymOp.cc:9
SymOp is the Coordinate representation of a symmetry operation it keeps fraction (FRAC) and Cartesian...
Definition: SymOp.hh:28
vector_type m_tau
Definition: SymOp.hh:162
EigenIndex Index
For long integer indexing:
const MasterSymGroup & master_group() const
const access of head group
Index ind_inverse() const
Get the operation index of the inverse of this operation, using the master_group's multiplication tab...
SymOp inverse() const
get the inverse of this SymOp
Definition: SymOp.cc:82
CASM::jsonParser & to_json_array(const Eigen::MatrixBase< Derived > &value, CASM::jsonParser &json)
Write Eigen Matrix with 1 row or 1 column to JSON array.
Definition: container.hh:191
double m_map_error
Definition: SymOp.hh:166
Eigen::Matrix3d matrix_type
Definition: SymOp.hh:31
jsonParser & put_obj()
Puts new empty JSON object.
Definition: jsonParser.hh:276
bool contains(const std::string &name) const
Return true if JSON object contains 'name'.
Definition: jsonParser.cc:500
void print(std::ostream &stream, const Eigen::Ref< const matrix_type > &c2fmat=matrix_type::Identity()) const
Definition: SymOp.cc:118
SymOp no_trans() const
Get copy of the SymOp without translation.
Definition: SymOp.cc:97
SymOp & apply_sym(const SymOp &op)
Definition: SymOp.cc:111
const vector_type & tau() const
Const access of the cartesian translation vector, 'tau'.
Definition: SymOp.hh:63
Index ind_prod(const SymOpRepresentation &RHS) const
SymOp & operator+=(const Eigen::Ref< const vector_type > &RHS)
Cartesian translation of the origin of the symmetry operation by Coordinate RHS.
Definition: SymOp.cc:66
bool almost_equal(const GenericCluster< CoordType > &LHS, const GenericCluster< CoordType > &RHS, double tol)
bool valid_index(Index i)